DAPartDeriv.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  Description:
7  Compute partial derivatives using the finite-difference method
8  with coloring
9 
10 \*---------------------------------------------------------------------------*/
11 
12 #ifndef DAPartDeriv_H
13 #define DAPartDeriv_H
14 
15 #include "runTimeSelectionTables.H"
16 #include "fvOptions.H"
17 #include "DAUtility.H"
18 #include "DAOption.H"
19 #include "DAIndex.H"
20 #include "DAModel.H"
21 #include "DAStateInfo.H"
22 #include "syncTools.H"
23 #include "DAObjFunc.H"
24 #include "DAJacCon.H"
25 #include "DAResidual.H"
26 
27 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
28 
29 namespace Foam
30 {
31 
32 /*---------------------------------------------------------------------------*\
33  Class DAPartDeriv Declaration
34 \*---------------------------------------------------------------------------*/
35 
37 {
38 
39 private:
41  DAPartDeriv(const DAPartDeriv&);
42 
44  void operator=(const DAPartDeriv&);
45 
46 protected:
48  const word modelType_;
49 
51  const fvMesh& mesh_;
52 
55 
57  const DAModel& daModel_;
58 
60  const DAIndex& daIndex_;
61 
64 
67 
69  const dictionary& allOptions_;
70 
72  HashTable<wordList> stateInfo_;
73 
75  void perturbStates(
76  const Vec jacConColors,
77  const Vec normStatePerturbVec,
78  const label colorI,
79  const scalar delta,
80  Vec wVec);
81 
83  void perturbBC(
84  const dictionary options,
85  const scalar delta);
86 
88  void setPartDerivMat(
89  const Vec resVec,
90  const Vec coloredColumn,
91  const label transposed,
92  Mat jacMat,
93  const scalar jacLowerBound=1e-30) const;
94 
96  void perturbAOA(
97  const dictionary options,
98  const scalar delta);
99 
102 
103 public:
105  TypeName("DAPartDeriv");
106 
107  // Declare run-time constructor selection table
109  autoPtr,
110  DAPartDeriv,
111  dictionary,
112  (const word modelType,
113  const fvMesh& mesh,
114  const DAOption& daOption,
115  const DAModel& daModel,
116  const DAIndex& daIndex,
117  const DAJacCon& daJacCon,
118  const DAResidual& daResidual),
119  (modelType,
120  mesh,
121  daOption,
122  daModel,
123  daIndex,
124  daJacCon,
125  daResidual));
126 
127  // Constructors
128 
129  //- Construct from components
130  DAPartDeriv(
131  const word modelType,
132  const fvMesh& mesh,
133  const DAOption& daOption,
134  const DAModel& daModel,
135  const DAIndex& daIndex,
136  const DAJacCon& daJacCon,
137  const DAResidual& daResidual);
138 
139  // Selectors
140 
141  //- Return a reference to the selected model
142  static autoPtr<DAPartDeriv> New(
143  const word modelType,
144  const fvMesh& mesh,
145  const DAOption& daOption,
146  const DAModel& daModel,
147  const DAIndex& daIndex,
148  const DAJacCon& daJacCon,
149  const DAResidual& daResidual);
150 
151  //- Destructor
152  virtual ~DAPartDeriv()
153  {
154  }
155 
156  // Member functions
157 
159  void clear()
160  {
161  MatDestroy(&dXvdFFDMat_);
162  }
163 
165  virtual void initializePartDerivMat(
166  const dictionary& options,
167  Mat jacMat) = 0;
168 
170  virtual void calcPartDerivMat(
171  const dictionary& options,
172  const Vec xvVec,
173  const Vec wVec,
174  Mat jacMat) = 0;
175 
177  void setdXvdFFDMat(const Mat dXvdFFDMat);
178 
180  void setNormStatePerturbVec(Vec* normStatePerturbVec);
181 };
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 } // End namespace Foam
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
Foam::DAJacCon
Definition: DAJacCon.H:34
DAOption.H
Foam::DAPartDeriv::clear
void clear()
clear members in parent and child objects
Definition: DAPartDeriv.H:159
DAIndex.H
Foam::DAOption
Definition: DAOption.H:29
DAUtility.H
Foam::DAPartDeriv::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAPartDeriv.H:60
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAPartDeriv::setPartDerivMat
void setPartDerivMat(const Vec resVec, const Vec coloredColumn, const label transposed, Mat jacMat, const scalar jacLowerBound=1e-30) const
set values for the partial derivative matrix
Definition: DAPartDeriv.C:165
Foam::DAPartDeriv::daJacCon_
const DAJacCon & daJacCon_
DAJacCon object.
Definition: DAPartDeriv.H:63
Foam::DAPartDeriv::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAPartDeriv.H:54
DAJacCon.H
Foam::DAPartDeriv::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAPartDeriv.H:51
DAResidual.H
DAModel.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAPartDeriv::stateInfo_
HashTable< wordList > stateInfo_
the stateInfo_ list from DAStateInfo object
Definition: DAPartDeriv.H:72
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAPartDeriv::New
static autoPtr< DAPartDeriv > New(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAJacCon &daJacCon, const DAResidual &daResidual)
Definition: DAPartDeriv.C:47
Foam::DAPartDeriv::~DAPartDeriv
virtual ~DAPartDeriv()
Definition: DAPartDeriv.H:152
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAPartDeriv::initializePartDerivMat
virtual void initializePartDerivMat(const dictionary &options, Mat jacMat)=0
initialize partial derivative matrix
Foam::DAPartDeriv::perturbAOA
void perturbAOA(const dictionary options, const scalar delta)
perturb the angle of attack
Definition: DAPartDeriv.C:375
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DAPartDeriv::TypeName
TypeName("DAPartDeriv")
Runtime type information.
Foam::DAPartDeriv::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, DAPartDeriv, dictionary,(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAJacCon &daJacCon, const DAResidual &daResidual),(modelType, mesh, daOption, daModel, daIndex, daJacCon, daResidual))
Foam::DAPartDeriv::perturbBC
void perturbBC(const dictionary options, const scalar delta)
perturb the values in the boundary condition
Definition: DAPartDeriv.C:266
DAStateInfo.H
Foam::DAPartDeriv::modelType_
const word modelType_
the name of the jacCon matrix
Definition: DAPartDeriv.H:48
daModel
DAModel daModel(mesh, daOption)
Foam::DAPartDeriv::setdXvdFFDMat
void setdXvdFFDMat(const Mat dXvdFFDMat)
setup dXvdFFD matrix
Definition: DAPartDeriv.C:498
Foam::DAPartDeriv
Definition: DAPartDeriv.H:36
Foam::DAPartDeriv::daResidual_
const DAResidual & daResidual_
DAResidual object.
Definition: DAPartDeriv.H:66
Foam::DAPartDeriv::perturbStates
void perturbStates(const Vec jacConColors, const Vec normStatePerturbVec, const label colorI, const scalar delta, Vec wVec)
perturb state variables given a color index
Definition: DAPartDeriv.C:98
Foam::DAPartDeriv::dXvdFFDMat_
Mat dXvdFFDMat_
volume mesh coordinates wrt the ffd point coordinate partials
Definition: DAPartDeriv.H:101
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAPartDeriv::setNormStatePerturbVec
void setNormStatePerturbVec(Vec *normStatePerturbVec)
setup the state normalization vector
Definition: DAPartDeriv.C:511
Foam::DAPartDeriv::calcPartDerivMat
virtual void calcPartDerivMat(const dictionary &options, const Vec xvVec, const Vec wVec, Mat jacMat)=0
compute the partial derivative matrix
DAObjFunc.H
Foam::DAPartDeriv::allOptions_
const dictionary & allOptions_
all the DAFoam option
Definition: DAPartDeriv.H:69
Foam::DAPartDeriv::daModel_
const DAModel & daModel_
DAModel object.
Definition: DAPartDeriv.H:57