DAModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
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 
46 #ifndef SolidDASolver
47 #include "DATurbulenceModel.H"
48 #include "DARadiationModel.H"
49 #endif
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class DAModel Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 class DAModel
60 {
61 
62 private:
64  DAModel(const DAModel&);
65 
67  void operator=(const DAModel&);
68 
70  const fvMesh& mesh_;
71 
73  const DAOption& daOption_;
74 
76  label hasTurbulenceModel_ = 0;
77 
79  label hasRadiationModel_ = 0;
80 
81 public:
83  DAModel(
84  const fvMesh& mesh,
85  const DAOption& daOption);
86 
88  virtual ~DAModel();
89 
91  void correctModelStates(wordList& modelStates) const;
92 
94  void correctStateResidualModelCon(List<List<word>>& stateCon) const;
95 
97  void addModelResidualCon(HashTable<List<List<word>>>& allCon) const;
98 
100  void calcResiduals(const dictionary& options);
101 
104 
107 
109  void getTurbProdTerm(scalarList& prodTerm) const;
110 
111 #ifndef SolidDASolver
114 #endif
115 
116 #ifdef CompressibleFlow
117  const fluidThermo& getThermo() const;
119 #endif
120 
121 };
122 
123 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124 
125 } // End namespace Foam
126 
127 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
128 
129 #endif
130 
131 // ************************************************************************* //
Foam::DAModel::getTurbProdTerm
void getTurbProdTerm(scalarList &prodTerm) const
return the value of the production term from the turbulence model
Definition: DAModel.C:273
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:73
Foam::DAModel::calcResiduals
void calcResiduals(const dictionary &options)
calculate the residuals for model state variables
Definition: DAModel.C:200
Foam::DAOption
Definition: DAOption.H:29
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAModel::~DAModel
virtual ~DAModel()
Destructor.
Definition: DAModel.C:29
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:127
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAModel::correctBoundaryConditions
void correctBoundaryConditions()
correct boundary conditions for model states
Definition: DAModel.C:227
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
DARadiationModel.H
Foam::DAModel::getDATurbulenceModel
const DATurbulenceModel & getDATurbulenceModel() const
get a reference to DATurbulenceModel
Definition: DAModel.C:185
Foam::DAModel::correctModelStates
void correctModelStates(wordList &modelStates) const
update the name in modelStates based on the selected physical model at runtime
Definition: DAModel.C:33
DATurbulenceModel.H
Foam::DAModel::updateIntermediateVariables
void updateIntermediateVariables()
update intermediate variables that are dependent on the model states
Definition: DAModel.C:250