DAPartDerivdFdAOA.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAPartDerivdFdAOA.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAPartDerivdFdAOA, 0);
16 addToRunTimeSelectionTable(DAPartDeriv, DAPartDerivdFdAOA, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  const word modelType,
21  const fvMesh& mesh,
22  const DAOption& daOption,
23  const DAModel& daModel,
24  const DAIndex& daIndex,
25  const DAJacCon& daJacCon,
26  const DAResidual& daResidual)
27  : DAPartDeriv(modelType,
28  mesh,
29  daOption,
30  daModel,
31  daIndex,
32  daJacCon,
33  daResidual)
34 {
35 }
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 
41  const dictionary& options,
42  Mat jacMat)
43 {
44  /*
45  Description:
46  Initialize jacMat
47 
48  Input:
49  options. This is not used
50  */
51 
52  // create dFdAOA
53  //MatCreate(PETSC_COMM_WORLD, jacMat);
54  MatSetSizes(
55  jacMat,
56  PETSC_DECIDE,
57  PETSC_DECIDE,
58  1,
59  1);
60  MatSetFromOptions(jacMat);
61  MatMPIAIJSetPreallocation(jacMat, 1, NULL, 1, NULL);
62  MatSeqAIJSetPreallocation(jacMat, 1, NULL);
63  //MatSetOption(jacMat, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
64  MatSetUp(jacMat);
65  MatZeroEntries(jacMat);
66  Info << "Partial derivative matrix created. " << mesh_.time().elapsedClockTime() << " s" << endl;
67 }
68 
70  const dictionary& options,
71  const Vec xvVec,
72  const Vec wVec,
73  Mat jacMat)
74 {
75  /*
76  Description:
77  Compute jacMat. Note for dFdAOA, we have only one column so we can do brute
78  force finite-difference there is no need to do coloring
79 
80  Input:
81  options.objFuncSubDictPart: the objFunc subDict, obtained from DAOption
82 
83  options.objFuncName: the name of the objective
84 
85  options.objFuncPart: the part of the objective
86 
87  options.varName: the name of the variable to perturb
88 
89  options.patchName: the name of the boundary patches to perturb, do not support
90  multiple patches yet
91 
92  options.flowAxis: the streamwise axis, aoa will be atan(U_normal/U_flow)
93 
94  options.normalAxis: the flow normal axis, aoa will be atan(U_normal/U_flow)
95 
96  xvVec: the volume mesh coordinate vector
97 
98  wVec: the state variable vector
99 
100  Output:
101  jacMat: the partial derivative matrix dFdAOA to compute
102  */
103 
104  word objFuncName, objFuncPart;
105  dictionary objFuncSubDictPart = options.subDict("objFuncSubDictPart");
106  options.readEntry<word>("objFuncName", objFuncName);
107  options.readEntry<word>("objFuncPart", objFuncPart);
108 
109  autoPtr<DAObjFunc> daObjFunc(
111  mesh_,
112  daOption_,
113  daModel_,
114  daIndex_,
115  daResidual_,
116  objFuncName,
117  objFuncPart,
118  objFuncSubDictPart)
119  .ptr());
120 
121  // zero all the matrices
122  MatZeroEntries(jacMat);
123 
124  dictionary mOptions;
125  mOptions.set("updateState", 1);
126  mOptions.set("updateMesh", 0);
127  scalar fRef = daObjFunc->masterFunction(mOptions, xvVec, wVec);
128 
129  scalar delta = daOption_.getSubDictOption<scalar>("adjPartDerivFDStep", "AOA");
130  scalar rDelta = 1.0 / delta;
131 
132  // perturb AOA
133  this->perturbAOA(options, delta);
134 
135  // compute object
136  scalar fNew = daObjFunc->masterFunction(mOptions, xvVec, wVec);
137 
138  scalar partDeriv = (fNew - fRef) * rDelta;
139  PetscScalar partDerivValue = 0.0;
140  assignValueCheckAD(partDerivValue, partDeriv);
141 
142  MatSetValue(jacMat, 0, 0, partDerivValue, INSERT_VALUES);
143 
144  // reset perturbation
145  this->perturbAOA(options, -1.0 * delta);
146  // call masterFunction again to reset the wVec to OpenFOAM field
147  fRef = daObjFunc->masterFunction(mOptions, xvVec, wVec);
148 
149  if (daOption_.getOption<label>("debug"))
150  {
151  Info << objFuncName << ": " << fRef << endl;
152  }
153 
154  MatAssemblyBegin(jacMat, MAT_FINAL_ASSEMBLY);
155  MatAssemblyEnd(jacMat, MAT_FINAL_ASSEMBLY);
156 }
157 
158 } // End namespace Foam
159 
160 // ************************************************************************* //
Foam::DAJacCon
Definition: DAJacCon.H:34
Foam::DAOption
Definition: DAOption.H:29
Foam::DAPartDerivdFdAOA::calcPartDerivMat
virtual void calcPartDerivMat(const dictionary &options, const Vec xvVec, const Vec wVec, Mat jacMat)
compute the partial derivative matrix
Definition: DAPartDerivdFdAOA.C:69
Foam::DAPartDeriv::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAPartDeriv.H:60
daOption
DAOption daOption(mesh, pyOptions_)
DAPartDerivdFdAOA.H
Foam::DAOption::getOption
classType getOption(const word key) const
get an option from subDict and key
Definition: DAOption.H:92
Foam::DAPartDeriv::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAPartDeriv.H:54
Foam::DAPartDerivdFdAOA::DAPartDerivdFdAOA
DAPartDerivdFdAOA(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAJacCon &daJacCon, const DAResidual &daResidual)
Definition: DAPartDerivdFdAOA.C:19
Foam::DAPartDeriv::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAPartDeriv.H:51
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAOption::getSubDictOption
classType getSubDictOption(const word subDict, const word subDictKey) const
get an dictionary option from subDict and key
Definition: DAOption.H:165
Foam::DAIndex
Definition: DAIndex.H:32
assignValueCheckAD
#define assignValueCheckAD(valA, valB)
Definition: DAMacroFunctions.H:103
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAPartDeriv::perturbAOA
void perturbAOA(const dictionary options, const scalar delta)
perturb the angle of attack
Definition: DAPartDeriv.C:375
Foam::DAPartDerivdFdAOA::initializePartDerivMat
virtual void initializePartDerivMat(const dictionary &options, Mat jacMat)
initialize partial derivative matrix
Definition: DAPartDerivdFdAOA.C:40
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::DAResidual
Definition: DAResidual.H:35
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
daModel
DAModel daModel(mesh, daOption)
Foam::DAPartDeriv
Definition: DAPartDeriv.H:36
Foam::DAPartDeriv::daResidual_
const DAResidual & daResidual_
DAResidual object.
Definition: DAPartDeriv.H:66
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAPartDeriv::daModel_
const DAModel & daModel_
DAModel object.
Definition: DAPartDeriv.H:57