DAModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAModel.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 // Constructors
16 DAModel::DAModel(
17  const fvMesh& mesh,
18  const DAOption& daOption)
19  : mesh_(mesh),
20  daOption_(daOption)
21 {
22  // check whether we have registered any physical models
23  hasTurbulenceModel_ = mesh.thisDb().foundObject<DATurbulenceModel>("DATurbulenceModel");
24  hasRadiationModel_ = mesh.thisDb().foundObject<DARadiationModel>("DARadiationModel");
25 }
26 
28 {
29 }
30 
31 void DAModel::correctModelStates(wordList& modelStates) const
32 {
33  /*
34  Description:
35  Update the name in modelStates based on the selected physical model at runtime
36 
37  Example:
38  In DAStateInfo, if the modelStates reads:
39 
40  modelStates = {"nut"}
41 
42  then for the SA model, calling correctModelStates(modelStates) will give:
43 
44  modelStates={"nuTilda"}
45 
46  while calling correctModelStates(modelStates) for the SST model will give
47 
48  modelStates={"k","omega"}
49 
50  We don't udpate the names for the radiation model becasue users are
51  supposed to set modelStates={"G"}
52  */
53 
54  // correct turbulence
55  if (hasTurbulenceModel_)
56  {
57  const DATurbulenceModel& daTurb =
58  mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel");
59  daTurb.correctModelStates(modelStates);
60  }
61 
62  // correct radiation
63  if (hasRadiationModel_)
64  {
65  // correct nothing because we should have register G for modelStates
66  }
67 }
68 
69 void DAModel::correctStateResidualModelCon(List<List<word>>& stateCon) const
70 {
71  /*
72  Description:
73  Update the original variable connectivity for the adjoint state
74  residuals in stateCon. Basically, we modify/add state variables based on the
75  original model variables defined in stateCon.
76 
77  Input:
78 
79  stateResCon: the connectivity levels for a state residual, defined in Foam::DAJacCon
80 
81  Example:
82  If stateCon reads:
83  stateCon=
84  {
85  {"U", "p", "nut"},
86  {"p"}
87  }
88 
89  For the SA turbulence model, calling this function for will get a new stateCon
90  stateCon=
91  {
92  {"U", "p", "nuTilda"},
93  {"p"}
94  }
95 
96  For the SST turbulence model, calling this function will give
97  stateCon=
98  {
99  {"U", "p", "k", "omega"},
100  {"p", "U"}
101  }
102  ***NOTE***: we add a extra level of U connectivity because nut is
103  related to grad(U), k, and omega in SST!
104  */
105 
106  // correct turbulence model states
107  if (hasTurbulenceModel_)
108  {
109  const DATurbulenceModel& daTurb =
110  mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel");
111  daTurb.correctStateResidualModelCon(stateCon);
112  }
113 
114  // correct radiation model states
115  if (hasRadiationModel_)
116  {
117  // correct nothing because we should have register G for modelStates
118  }
119 }
120 
121 void DAModel::addModelResidualCon(HashTable<List<List<word>>>& allCon) const
122 {
123  /*
124  Description:
125  Add the connectivity levels for all physical model residuals to allCon
126 
127  Input:
128  allCon: the connectivity levels for all state residual, defined in DAJacCon
129 
130  Example:
131  If stateCon reads:
132  allCon=
133  {
134  "URes":
135  {
136  {"U", "p", "nut"},
137  {"p"}
138  }
139  }
140 
141  For the SA turbulence model, calling this function for will get a new stateCon,
142  something like this:
143  allCon=
144  {
145  "URes":
146  {
147  {"U", "p", "nuTilda"},
148  {"p"}
149  },
150  "nuTildaRes":
151  {
152  {"U", "phi", "nuTilda"},
153  {"U"}
154  }
155  }
156 
157  */
158 
159  // add turbulence model state residuals
160  if (hasTurbulenceModel_)
161  {
162  const DATurbulenceModel& daTurb =
163  mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel");
164  daTurb.addModelResidualCon(allCon);
165  }
166 
167  // add radiation model state residuals
168  if (hasRadiationModel_)
169  {
170  const DARadiationModel& daRadiation =
171  mesh_.thisDb().lookupObject<DARadiationModel>("DARadiationModel");
172  daRadiation.addModelResidualCon(allCon);
173  }
174 }
175 
177 {
178  /*
179  Description:
180  Return the Foam::DATurbulence object
181  */
182  if (!hasTurbulenceModel_)
183  {
184  FatalErrorIn("getDATurbulenceModel") << "DATurbulenceModel not found in mesh.thisDb()! "
185  << abort(FatalError);
186  }
187  return mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel");
188 }
189 
190 void DAModel::calcResiduals(const dictionary& options)
191 {
192  /*
193  Description:
194  Calculate the residuals for model state variables
195 
196  Input:
197  options.isPC: whether to compute simplfied residual for preconditoner
198  See the child classes in Foam::DATurbulenceModel for details
199  */
200 
201  if (hasTurbulenceModel_)
202  {
203  DATurbulenceModel& daTurb = const_cast<DATurbulenceModel&>(
204  mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel"));
205  daTurb.calcResiduals(options);
206  }
207 
208  if (hasRadiationModel_)
209  {
210  // not implemented
211  }
212 }
213 
215 {
216  /*
217  Description:
218  Update the boundary conditions for model variables. Check the child classes
219  in Foam::DATurbulenceModel for details
220  */
221 
222  if (hasTurbulenceModel_)
223  {
224  DATurbulenceModel& daTurb = const_cast<DATurbulenceModel&>(
225  mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel"));
226  daTurb.correctBoundaryConditions();
227  }
228 
229  if (hasRadiationModel_)
230  {
231  }
232 }
233 
235 {
236  /*
237  Description:
238  Update the intermediate variables for models, e.g., update nut based on nuTilda
239  Check the child classes in Foam::DATurbulenceModel for details
240  */
241 
242  if (hasTurbulenceModel_)
243  {
244  DATurbulenceModel& daTurb = const_cast<DATurbulenceModel&>(
245  mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel"));
247  }
248 
249  if (hasRadiationModel_)
250  {
251  }
252 }
253 
254 void DAModel::getTurbProdTerm(volScalarField& prodTerm) const
255 {
256  /*
257  Description:
258  Return the value of the production term from the turbulence model
259  */
260 
261  if (hasTurbulenceModel_)
262  {
263  DATurbulenceModel& daTurb = const_cast<DATurbulenceModel&>(
264  mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel"));
265  daTurb.getTurbProdTerm(prodTerm);
266  }
267 
268  if (hasRadiationModel_)
269  {
270  }
271 }
272 
273 void DAModel::getTurbProdOverDestruct(volScalarField& PoD) const
274 {
275  /*
276  Description:
277  Return the value of the production/destruction term from the turbulence model
278  */
279 
280  if (hasTurbulenceModel_)
281  {
282  DATurbulenceModel& daTurb = const_cast<DATurbulenceModel&>(
283  mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel"));
284  daTurb.getTurbProdOverDestruct(PoD);
285  }
286 
287  if (hasRadiationModel_)
288  {
289  }
290 }
291 
292 void DAModel::getTurbConvOverProd(volScalarField& CoP) const
293 {
294  /*
295  Description:
296  return the value of the convective over production term from the turbulence model
297  */
298 
299  if (hasTurbulenceModel_)
300  {
301  DATurbulenceModel& daTurb = const_cast<DATurbulenceModel&>(
302  mesh_.thisDb().lookupObject<DATurbulenceModel>("DATurbulenceModel"));
303  daTurb.getTurbConvOverProd(CoP);
304  }
305 
306  if (hasRadiationModel_)
307  {
308  }
309 }
310 
311 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
312 
313 } // End namespace Foam
314 
315 // ************************************************************************* //
Foam::DARadiationModel::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const =0
Foam::DATurbulenceModel::getTurbProdTerm
virtual void getTurbProdTerm(volScalarField &prodTerm) const
return the value of the production term from the turbulence model
Definition: DATurbulenceModel.C:540
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::DATurbulenceModel::updateIntermediateVariables
virtual void updateIntermediateVariables()=0
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Foam::DAOption
Definition: DAOption.H:29
Foam::DATurbulenceModel::correctModelStates
virtual void correctModelStates(wordList &modelStates) const =0
update the turbulence state for DAStateInfo::regStates_
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
DAModel.H
Foam::DATurbulenceModel::correctBoundaryConditions
virtual void correctBoundaryConditions()=0
update turbulence variable boundary values
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::DATurbulenceModel::calcResiduals
virtual void calcResiduals(const dictionary &options)=0
compute the turbulence residuals
Foam
Definition: checkGeometry.C:32
Foam::DATurbulenceModel::getTurbConvOverProd
virtual void getTurbConvOverProd(volScalarField &CoP) const
return the value of the convective over production term from the turbulence model
Definition: DATurbulenceModel.C:564
Foam::DAModel::correctBoundaryConditions
void correctBoundaryConditions()
correct boundary conditions for model states
Definition: DAModel.C:214
Foam::DATurbulenceModel::addModelResidualCon
virtual void addModelResidualCon(HashTable< List< List< word >>> &allCon) const =0
add the model residual connectivity to stateCon
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
Foam::DATurbulenceModel::correctStateResidualModelCon
virtual void correctStateResidualModelCon(List< List< word >> &stateCon) const =0
update the original variable connectivity for the adjoint state residuals in stateCon
Foam::DARadiationModel
Definition: DARadiationModel.H:27
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
Foam::DATurbulenceModel::getTurbProdOverDestruct
virtual void getTurbProdOverDestruct(volScalarField &PoD) const
return the ratio of the production over destruction term from the turbulence model
Definition: DATurbulenceModel.C:552
Foam::DAModel::updateIntermediateVariables
void updateIntermediateVariables()
update intermediate variables that are dependent on the model states
Definition: DAModel.C:234