DAResidual.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAResidual.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
16 
17 defineTypeNameAndDebug(DAResidual, 0);
18 defineRunTimeSelectionTable(DAResidual, dictionary);
19 
20 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
21 
22 DAResidual::DAResidual(
23  const word modelType,
24  const fvMesh& mesh,
25  const DAOption& daOption,
26  const DAModel& daModel,
27  const DAIndex& daIndex)
28  : mesh_(mesh),
29  daOption_(daOption),
30  daModel_(daModel),
31  daIndex_(daIndex),
32  daField_(mesh, daOption, daModel, daIndex)
33 {
34 }
35 
36 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
37 
38 autoPtr<DAResidual> DAResidual::New(
39  const word modelType,
40  const fvMesh& mesh,
41  const DAOption& daOption,
42  const DAModel& daModel,
43  const DAIndex& daIndex)
44 {
45  // standard setup for runtime selectable classes
46 
47  if (daOption.getAllOptions().lookupOrDefault<label>("debug", 0))
48  {
49  Info << "Selecting " << modelType << " for DAResidual" << endl;
50  }
51 
52  dictionaryConstructorTable::iterator cstrIter =
53  dictionaryConstructorTablePtr_->find(modelType);
54 
55  // if the solver name is not found in any child class, print an error
56  if (cstrIter == dictionaryConstructorTablePtr_->end())
57  {
58  FatalErrorIn(
59  "DAResidual::New"
60  "("
61  " const word,"
62  " const fvMesh&,"
63  " const DAOption&,"
64  " const DAModel&,"
65  " const DAIndex&"
66  ")")
67  << "Unknown DAResidual type "
68  << modelType << nl << nl
69  << "Valid DAResidual types:" << endl
70  << dictionaryConstructorTablePtr_->sortedToc()
71  << exit(FatalError);
72  }
73 
74  // child class found
75  return autoPtr<DAResidual>(
76  cstrIter()(modelType, mesh, daOption, daModel, daIndex));
77 }
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
82  const dictionary& options,
83  const Vec xvVec,
84  const Vec wVec,
85  Vec resVec)
86 {
87  /*
88  Description:
89  A master function that takes the volume mesh points and state variable vecs
90  as input, and compute the residual vector.
91 
92  Input:
93  options.updateState: whether to assign the values in wVec to the state
94  variables of the OpenFOAM fields (e.g., U, p). This will also update boundary conditions
95  and update all intermediate variables that are dependent on the state
96  variables.
97 
98  options.updateMesh: whether to assign the values in xvVec to the OpenFOAM mesh
99  coordinates in Foam::fvMesh. This will also call mesh.movePoints() to update
100  all the mesh metrics such as mesh volume, cell centers, mesh surface area, etc.
101 
102  options.setResVec: whether to assign the residuals (e.g., URes_) from the OpenFOAM fields
103  to the resVec. If set to 0, nothing will be written to resVec
104 
105  options.isPC: whether to compute residual for constructing PC matrix
106 
107  xvVec: the volume coordinates vector (flatten)
108 
109  wVec: the state variable vector
110 
111  Output:
112  resVec: the residual vector
113 
114  NOTE1: the ordering of the xvVec, wVec, and resVec depends on adjStateOrdering, see
115  Foam::DAIndex for details
116 
117  NOTE2: the calcResiduals function will be implemented in the child classes
118  */
119 
120  VecZeroEntries(resVec);
121 
122  DAModel& daModel = const_cast<DAModel&>(daModel_);
123 
124  label updateState = options.getLabel("updateState");
125 
126  label updateMesh = options.getLabel("updateMesh");
127 
128  label setResVec = options.getLabel("setResVec");
129 
130  if (updateMesh)
131  {
132  daField_.pointVec2OFMesh(xvVec);
133  }
134 
135  if (updateState)
136  {
138 
139  // now update intermediate states and boundry conditions
142  daModel.correctBoundaryConditions();
143  daModel.updateIntermediateVariables();
144  // if there are special boundary conditions, apply special treatment
146  }
147 
148  this->calcResiduals(options);
149  daModel.calcResiduals(options);
150 
151  if (setResVec)
152  {
153  // asssign the openfoam residual field to resVec
154  daField_.ofResField2ResVec(resVec);
155  }
156 }
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace Foam
161 
162 // ************************************************************************* //
Foam::DAResidual::correctBoundaryConditions
virtual void correctBoundaryConditions()=0
update the boundary condition for all the states in the selected solver
Foam::DAOption
Definition: DAOption.H:29
Foam::DAResidual::updateIntermediateVariables
virtual void updateIntermediateVariables()=0
update any intermdiate variables that are dependent on state variables and are used in calcResiduals
daOption
DAOption daOption(mesh, pyOptions_)
DAResidual.H
Foam::DAResidual::daField_
DAField daField_
DAField object.
Definition: DAResidual.H:59
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(DAFvSource, dictionary)
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAResidual::New
static autoPtr< DAResidual > New(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidual.C:38
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::DAField::stateVec2OFField
void stateVec2OFField(const Vec stateVec) const
assign the fields in OpenFOAM based on the state vector
Definition: DAField.C:375
Foam::DAResidual::daModel_
const DAModel & daModel_
DAModel object.
Definition: DAResidual.H:53
Foam::DAResidual::masterFunction
void masterFunction(const dictionary &options, const Vec xvVec, const Vec wVec, Vec resVec)
the master function that compute the residual vector given the state and point vectors
Definition: DAResidual.C:81
Foam::DAField::ofResField2ResVec
void ofResField2ResVec(Vec resVec) const
assign the residual vector based on the residual field in OpenFOAM
Definition: DAField.C:591
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
daModel
DAModel daModel(mesh, daOption)
Foam::DAField::pointVec2OFMesh
void pointVec2OFMesh(const Vec xvVec) const
assign the points in fvMesh of OpenFOAM based on the point vector
Definition: DAField.C:506
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAField::specialBCTreatment
void specialBCTreatment()
apply special treatment for boundary conditions
Definition: DAField.C:812
Foam::DAResidual::calcResiduals
virtual void calcResiduals(const dictionary &options)=0
compute residuals