DAObjFunc.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  Description:
7  Compute object functions and their derivatives
8 
9 \*---------------------------------------------------------------------------*/
10 
11 #ifndef DAObjFunc_H
12 #define DAObjFunc_H
13 
14 #include "runTimeSelectionTables.H"
15 #include "fvOptions.H"
16 #include "DAOption.H"
17 #include "DAModel.H"
18 #include "DAIndex.H"
19 #include "topoSetSource.H"
20 #include "topoSet.H"
21 #include "DAField.H"
22 #include "DAResidual.H"
23 
24 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
25 
26 namespace Foam
27 {
28 
29 /*---------------------------------------------------------------------------*\
30  Class DAObjFunc Declaration
31 \*---------------------------------------------------------------------------*/
32 
33 class DAObjFunc
34 {
35 
36 private:
38  DAObjFunc(const DAObjFunc&);
39 
41  void operator=(const DAObjFunc&);
42 
43 protected:
45  const fvMesh& mesh_;
46 
49 
51  const DAModel& daModel_;
52 
54  const DAIndex& daIndex_;
55 
58 
61 
64 
67 
69  const dictionary& objFuncDict_;
70 
73 
76 
79 
81  scalarList objFuncFaceValues_;
82 
84  scalarList objFuncCellValues_;
85 
87  scalar scale_;
88 
90  scalar objFuncValue_;
91 
93  List<List<word>> objFuncConInfo_;
94 
95 public:
97  TypeName("DAObjFunc");
98 
99  // Declare run-time constructor selection table
101  autoPtr,
102  DAObjFunc,
103  dictionary,
104  (
105  const fvMesh& mesh,
106  const DAOption& daOption,
107  const DAModel& daModel,
108  const DAIndex& daIndex,
109  const DAResidual& daResidual,
110  const word objFuncName,
111  const word objFuncPart,
112  const dictionary& objFuncDict),
113  (
114  mesh,
115  daOption,
116  daModel,
117  daIndex,
118  daResidual,
119  objFuncName,
120  objFuncPart,
121  objFuncDict));
122 
123  // Constructors
124 
125  //- Construct from components
126  DAObjFunc(
127  const fvMesh& mesh,
128  const DAOption& daOption,
129  const DAModel& daModel,
130  const DAIndex& daIndex,
131  const DAResidual& daResidual,
132  const word objFuncName,
133  const word objFuncPart,
134  const dictionary& objFuncDict);
135 
136  // Selectors
137 
138  //- Return a reference to the selected model
139  static autoPtr<DAObjFunc> New(
140  const fvMesh& mesh,
141  const DAOption& daOption,
142  const DAModel& daModel,
143  const DAIndex& daIndex,
144  const DAResidual& daResidual,
145  const word objFuncName,
146  const word objFuncPart,
147  const dictionary& objFuncDict);
148 
149  //- Destructor
150  virtual ~DAObjFunc()
151  {
152  }
153 
155  void clear()
156  {
157  objFuncFaceSources_.clear();
158  objFuncCellSources_.clear();
159  objFuncFaceValues_.clear();
160  objFuncCellValues_.clear();
161  objFuncConInfo_.clear();
162  }
163 
165  void calcObjFuncSources(
166  labelList& faceSources,
167  labelList& cellSources);
168 
170  virtual void calcObjFunc(
171  const labelList& objFuncFaceSources,
172  const labelList& objFuncCellSources,
173  scalarList& objFuncFaceValues,
174  scalarList& objFuncCellValues,
175  scalar& objFuncValue) = 0;
176 
178  scalar getObjFuncValue();
179 
181  scalar masterFunction(
182  const dictionary& options,
183  const Vec xvVec,
184  const Vec wVec);
185 
188  {
189  return objFuncName_;
190  }
191 
194  {
195  return objFuncPart_;
196  }
197 
200  {
201  return objFuncType_;
202  }
203 
205  const labelList& getObjFuncFaceSources() const
206  {
207  return objFuncFaceSources_;
208  }
209 
211  const labelList& getObjFuncCellSources() const
212  {
213  return objFuncCellSources_;
214  }
215 
217  const scalarList& getObjFuncFaceValues() const
218  {
219  return objFuncFaceValues_;
220  }
221 
223  const scalarList& getObjFuncCellValues() const
224  {
225  return objFuncCellValues_;
226  }
227 
229  const List<List<word>>& getObjFuncConInfo() const
230  {
231  return objFuncConInfo_;
232  }
233 
235  scalar expSumKS = 1.0;
236 
238  label calcRefCoeffs = 1;
239 };
240 
241 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242 
243 } // End namespace Foam
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 
247 #endif
248 
249 // ************************************************************************* //
DAOption.H
Foam::DAObjFunc::daResidual_
const DAResidual & daResidual_
DAResidual object.
Definition: DAObjFunc.H:57
Foam::DAObjFunc::objFuncType_
word objFuncType_
the type of the objective function
Definition: DAObjFunc.H:66
Foam::DAObjFunc::clear
void clear()
clear up members
Definition: DAObjFunc.H:155
Foam::DAObjFunc::getObjFuncCellValues
const scalarList & getObjFuncCellValues() const
return DAObjFunc::objFuncCellValues_
Definition: DAObjFunc.H:223
Foam::DAObjFunc::objFuncCellSources_
labelList objFuncCellSources_
a sorted list of all cell sources for the objective function
Definition: DAObjFunc.H:78
Foam::DAObjFunc::getObjFuncFaceValues
const scalarList & getObjFuncFaceValues() const
return DAObjFunc::objFuncFaceValues_
Definition: DAObjFunc.H:217
Foam::DAObjFunc::daModel_
const DAModel & daModel_
DAModel object.
Definition: DAObjFunc.H:51
Foam::DAObjFunc::objFuncName_
word objFuncName_
the name of the objective function
Definition: DAObjFunc.H:60
Foam::DAObjFunc::getObjFuncCellSources
const labelList & getObjFuncCellSources() const
return DAObjFunc::objFuncCellSources_
Definition: DAObjFunc.H:211
DAIndex.H
Foam::DAObjFunc::daField_
DAField daField_
DAField object.
Definition: DAObjFunc.H:72
Foam::DAOption
Definition: DAOption.H:29
Foam::DAObjFunc::getObjFuncFaceSources
const labelList & getObjFuncFaceSources() const
return DAObjFunc::objFuncFaceSources_
Definition: DAObjFunc.H:205
Foam::DAObjFunc::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAObjFunc.H:54
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAObjFunc::objFuncCellValues_
scalarList objFuncCellValues_
value of the obj function in the discrete cells
Definition: DAObjFunc.H:84
Foam::DAObjFunc::calcObjFuncSources
void calcObjFuncSources(labelList &faceSources, labelList &cellSources)
calculate DAObjFunc::objFuncFaceSources_ and DAObjFunc::objFuncCellSources_
Definition: DAObjFunc.C:123
Foam::DAObjFunc::calcObjFunc
virtual void calcObjFunc(const labelList &objFuncFaceSources, const labelList &objFuncCellSources, scalarList &objFuncFaceValues, scalarList &objFuncCellValues, scalar &objFuncValue)=0
calculate the value of objective function
Foam::DAObjFunc::getObjFuncPart
word getObjFuncPart()
return the part of objective function
Definition: DAObjFunc.H:193
DAResidual.H
Foam::DAObjFunc::getObjFuncValue
scalar getObjFuncValue()
calculate the value of objective function
Definition: DAObjFunc.C:290
Foam::DAObjFunc::getObjFuncType
word getObjFuncType()
return the part of objective function
Definition: DAObjFunc.H:199
DAModel.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAObjFunc::masterFunction
scalar masterFunction(const dictionary &options, const Vec xvVec, const Vec wVec)
the master function to compute objective function given the state and point vectors
Definition: DAObjFunc.C:228
Foam::DAObjFunc::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAObjFunc.H:45
Foam::DAObjFunc::getObjFuncConInfo
const List< List< word > > & getObjFuncConInfo() const
return DAObjFunc::objFuncConInfo_
Definition: DAObjFunc.H:229
Foam::DAObjFunc::~DAObjFunc
virtual ~DAObjFunc()
Definition: DAObjFunc.H:150
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAObjFunc::calcRefCoeffs
label calcRefCoeffs
whether to compute reference coefficients for special objFunc treatment such as totalPressureRatio
Definition: DAObjFunc.H:238
Foam::DAField
Definition: DAField.H:36
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAObjFunc::New
static autoPtr< DAObjFunc > New(const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAResidual &daResidual, const word objFuncName, const word objFuncPart, const dictionary &objFuncDict)
Definition: DAObjFunc.C:62
Foam::DAObjFunc::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAObjFunc.H:48
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DAObjFunc::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, DAObjFunc, dictionary,(const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAResidual &daResidual, const word objFuncName, const word objFuncPart, const dictionary &objFuncDict),(mesh, daOption, daModel, daIndex, daResidual, objFuncName, objFuncPart, objFuncDict))
Foam::DAObjFunc::TypeName
TypeName("DAObjFunc")
Runtime type information.
Foam::DAObjFunc::objFuncPart_
word objFuncPart_
the part of the objective function
Definition: DAObjFunc.H:63
Foam::DAObjFunc::scale_
scalar scale_
scale of the objective function
Definition: DAObjFunc.H:87
daModel
DAModel daModel(mesh, daOption)
Foam::DAObjFunc::objFuncConInfo_
List< List< word > > objFuncConInfo_
the connectivity information for the objective function
Definition: DAObjFunc.H:93
Foam::DAObjFunc::objFuncFaceValues_
scalarList objFuncFaceValues_
value of the obj function on the discrete faces
Definition: DAObjFunc.H:81
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAObjFunc::objFuncDict_
const dictionary & objFuncDict_
dictionary containing the information for the objective function
Definition: DAObjFunc.H:69
DAField.H
Foam::DAObjFunc::getObjFuncName
word getObjFuncName()
return the name of objective function
Definition: DAObjFunc.H:187
Foam::DAObjFunc::objFuncValue_
scalar objFuncValue_
the value of the objective function
Definition: DAObjFunc.H:90
Foam::DAObjFunc::objFuncFaceSources_
labelList objFuncFaceSources_
a sorted list of all face sources for the objective function
Definition: DAObjFunc.H:75
Foam::DAObjFunc
Definition: DAObjFunc.H:33
Foam::DAObjFunc::expSumKS
scalar expSumKS
expSumKS stores sum[exp(coeffKS*x_i)] for KS function which will be used to scale dFdW
Definition: DAObjFunc.H:235