DAModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  Description:
7  DAModel is a general class for handling physical models in DAFoam.
8  The constructor for DAModel does not explicitly contain objects of
9  physical models. Instead, it will look up any model class, e.g.,
10  DATurbulenceModel, DARaditionModel, that has been registred in fvMesh.
11  This enables the flexibility such that we can simultaneously have
12  multiple models for adjoint, e.g., turbulence model + radiation model.
13 
14  NOTE 1: in adjoint functions, we should call functions in DAModel to handle
15  general model operation, e.g., call daModel.correctBoundaryCondition(),
16  instead of calling daTurbulenceModel.correctBoundaryCondition(). Because
17  daModel.correctBoundaryCondition() will call the correctBoundaryCondition()
18  in daTurbulenceModel and daRadiationModel, if these two models are used.
19 
20  NOTE 2: Because of the above feature, we need to initialize any physical
21  model BEFORE initializing the DAModel, e.g.,
22 
23  // initialize DATurbulence and automatically register it to fvMesh
24  autoPtr<DATurbulenceModel> daTurb(DATurbulenceModel::New(mesh));
25 
26  ............ // codes to initialize any other physical models
27 
28  // now we can initialize DAModel and lookup DATurbulenceModel in mesh
29  DAModel daModel(mesh);
30 
31  Otherwise, the physical models will NOT be picked up by the DAModel
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef DAModel_H
36 #define DAModel_H
37 
38 #include "fvOptions.H"
39 #include "DAOption.H"
40 #include "fvc.H"
41 #include "fvm.H"
42 #include "surfaceFields.H"
43 #include "geometricOneField.H"
44 #include "zeroGradientFvPatchField.H"
45 #include "DATurbulenceModel.H"
46 #include "DARadiationModel.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 /*---------------------------------------------------------------------------*\
54  Class DAModel Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 class DAModel
58 {
59 
60 private:
62  DAModel(const DAModel&);
63 
65  void operator=(const DAModel&);
66 
68  const fvMesh& mesh_;
69 
71  const DAOption& daOption_;
72 
74  label hasTurbulenceModel_ = 0;
75 
77  label hasRadiationModel_ = 0;
78 
79 public:
81  DAModel(
82  const fvMesh& mesh,
83  const DAOption& daOption);
84 
86  virtual ~DAModel();
87 
89  void correctModelStates(wordList& modelStates) const;
90 
92  void correctStateResidualModelCon(List<List<word>>& stateCon) const;
93 
95  void addModelResidualCon(HashTable<List<List<word>>>& allCon) const;
96 
98  void calcResiduals(const dictionary& options);
99 
102 
105 
107  void getTurbProdTerm(volScalarField& prodTerm) const;
108 
110  void getTurbProdOverDestruct(volScalarField& PoD) const;
111 
113  void getTurbConvOverProd(volScalarField& CoP) const;
114 
117 };
118 
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120 
121 } // End namespace Foam
122 
123 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124 
125 #endif
126 
127 // ************************************************************************* //
DAOption.H
Foam::DAModel::correctStateResidualModelCon
void correctStateResidualModelCon(List< List< word >> &stateCon) const
update the original variable connectivity for the adjoint state residuals in stateCon
Definition: DAModel.C:69
Foam::DAModel::calcResiduals
void calcResiduals(const dictionary &options)
calculate the residuals for model state variables
Definition: DAModel.C:190
Foam::DAOption
Definition: DAOption.H:29
Foam::DAModel::getTurbProdTerm
void getTurbProdTerm(volScalarField &prodTerm) const
return the value of the production term from the turbulence model
Definition: DAModel.C:254
Foam::DAModel::~DAModel
virtual ~DAModel()
Destructor.
Definition: DAModel.C:27
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAModel::addModelResidualCon
void addModelResidualCon(HashTable< List< List< word >>> &allCon) const
add the model residual connectivity to stateCon
Definition: DAModel.C:121
Foam::DAModel
Definition: DAModel.H:57
Foam
Definition: checkGeometry.C:32
Foam::DAModel::correctBoundaryConditions
void correctBoundaryConditions()
correct boundary conditions for model states
Definition: DAModel.C:214
Foam::DAModel::getTurbConvOverProd
void getTurbConvOverProd(volScalarField &CoP) const
return the value of the convective over production term from the turbulence model
Definition: DAModel.C:292
Foam::DAModel::getTurbProdOverDestruct
void getTurbProdOverDestruct(volScalarField &PoD) const
return the value of the destruction term from the turbulence model
Definition: DAModel.C:273
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:52
DARadiationModel.H
Foam::DAModel::getDATurbulenceModel
const DATurbulenceModel & getDATurbulenceModel() const
get a reference to DATurbulenceModel
Definition: DAModel.C:176
Foam::DAModel::correctModelStates
void correctModelStates(wordList &modelStates) const
update the name in modelStates based on the selected physical model at runtime
Definition: DAModel.C:31
DATurbulenceModel.H
Foam::DAModel::updateIntermediateVariables
void updateIntermediateVariables()
update intermediate variables that are dependent on the model states
Definition: DAModel.C:234