DAStateInfoRhoSimpleCFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAStateInfoRhoSimpleCFoam, 0);
16 addToRunTimeSelectionTable(DAStateInfo, DAStateInfoRhoSimpleCFoam, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  const word modelType,
21  const fvMesh& mesh,
22  const DAOption& daOption,
23  const DAModel& daModel)
24  : DAStateInfo(modelType, mesh, daOption, daModel)
25 {
26  /*
27  Description:
28  Register the names of state variables
29  NOTE:
30  For model variables, such as turbulence model, register specific names
31  For example, register "nut" to modelStates for RANS turbulence models,
32  Then, we will call correctModelStates(stateInfo_["modelStates"]) to modify
33  "nut" based on the selected turbulence model. For example, for SA model,
34  correctModelStates will just replace "nut" with "nuTilda", for SST model,
35  it will replace "nut" with "k" and append "omega" to modelStates.
36  In other words, the model variables will be modified based on the selected
37  models at runtime.
38  */
39 
40  stateInfo_["volScalarStates"].append("p");
41  stateInfo_["volScalarStates"].append("T");
42  stateInfo_["modelStates"].append("nut");
43  stateInfo_["volVectorStates"].append("U");
44  stateInfo_["surfaceScalarStates"].append("phi");
45 
46  // correct the names for model states based on the selected physical model at runtime
47  daModel.correctModelStates(stateInfo_["modelStates"]);
48 
49  /*
50  Description:
51  Adjoint state connectivity info, numbers denote the level of connectivity
52  N/A means this state does not connect to the corrsponding residual
53 
54  U T p nut phi
55  URes 2 2 1 1 0
56  TRes 2 2 2 1 0
57  pRes 3 2 2 2 1
58  phiRes 2 2 1 1 0
59 
60  ******************************** NOTE 1 **********************************
61  One does not need to specify connectivity for each physical model, set the
62  connectivity for original variables instead. For example, for turbulence models,
63  set nut. Then, how is nut connected to the other turbulence states will be
64  set in the DAModel class. This is done by calling correctStateResidualModelCon.
65  For example, for SA model we just replace nut with nuTilda, for SST model, we need
66  to add extract connectivity since nut depends on grad(U), k, and omega. We need
67  to do this for other pyhsical models such as radiation models.
68  **************************************************************************
69 
70  ******************************** NOTE 2 **********************************
71  Do not specify physical model connectivity here, because they will be added
72  by calling addModelResidualCon. For example, for the SA turbulence
73  model, it will add the nuTildaRes to stateResConInfo_ and setup
74  its connectivity automatically.
75  **************************************************************************
76 
77  */
78 
79  stateResConInfo_.set(
80  "URes",
81  {
82  {"U", "p", "T", "nut", "phi"}, // lv0
83  {"U", "p", "T", "nut"}, // lv1
84  {"U", "T"} // lv2
85  });
86 
87  stateResConInfo_.set(
88  "TRes",
89  {
90  {"U", "p", "T", "nut", "phi"}, // lv0
91  {"U", "p", "T", "nut"}, // lv1
92  {"U", "p", "T"} // lv2
93  });
94 
95  stateResConInfo_.set(
96  "pRes",
97  {
98  {"U", "p", "T", "nut", "phi"}, // lv0
99  {"U", "p", "T", "nut", "phi"}, // lv1
100  {"U", "p", "T", "nut"}, // lv2
101  {"U"} // lv3
102  });
103 
104  stateResConInfo_.set(
105  "phiRes",
106  {
107  {"U", "p", "T", "nut", "phi"}, // lv0
108  {"U", "p", "T", "nut"}, // lv1
109  {"U", "T"}, // lv2
110  });
111 
112  // need to correct connectivity for physical models for each residual
113  daModel.correctStateResidualModelCon(stateResConInfo_["URes"]);
114  daModel.correctStateResidualModelCon(stateResConInfo_["TRes"]);
115  daModel.correctStateResidualModelCon(stateResConInfo_["pRes"]);
116  daModel.correctStateResidualModelCon(stateResConInfo_["phiRes"]);
117 
118  // add physical model residual connectivity
119  daModel.addModelResidualCon(stateResConInfo_);
120 
121 }
122 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123 
124 } // End namespace Foam
125 
126 // ************************************************************************* //
Foam::DAOption
Definition: DAOption.H:29
daOption
DAOption daOption(mesh, pyOptions_)
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
DAStateInfoRhoSimpleCFoam.H
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAStateInfoRhoSimpleCFoam::DAStateInfoRhoSimpleCFoam
DAStateInfoRhoSimpleCFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel)
Definition: DAStateInfoRhoSimpleCFoam.C:19
Foam::DAStateInfo::stateResConInfo_
HashTable< List< List< word > > > stateResConInfo_
table to specify how the states are connected to the residuals for a given solver
Definition: DAStateInfo.H:54
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
daModel
DAModel daModel(mesh, daOption)
Foam::DAStateInfo
Definition: DAStateInfo.H:30
Foam::DAStateInfo::stateInfo_
HashTable< wordList > stateInfo_
registered states 1st hash: solverName, 2nd hash: fieldType, 3nd list, stateNames
Definition: DAStateInfo.H:51