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