DAPartDerivdRdACTD.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAPartDerivdRdACTD.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAPartDerivdRdACTD, 0);
16 addToRunTimeSelectionTable(DAPartDeriv, DAPartDerivdRdACTD, 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 
49  // now initialize the memory for the jacobian itself
50  label localSize = daIndex_.nLocalAdjointStates;
51 
52  // create dRdACTDT
53  //MatCreate(PETSC_COMM_WORLD, jacMat);
54  MatSetSizes(
55  jacMat,
56  localSize,
57  PETSC_DECIDE,
58  PETSC_DETERMINE,
59  nActDVs_);
60  MatSetFromOptions(jacMat);
61  MatMPIAIJSetPreallocation(jacMat, nActDVs_, NULL, nActDVs_, NULL);
62  MatSeqAIJSetPreallocation(jacMat, nActDVs_, 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. We use brute-force finite-difference
78 
79  Input:
80 
81  options.isPC: whether to compute the jacMat for preconditioner
82 
83  xvVec: the volume mesh coordinate vector
84 
85  wVec: the state variable vector
86 
87  Output:
88  jacMat: the partial derivative matrix dRdACTD to compute
89  */
90 
91  word actuatorName = options.getWord("actuatorName");
92 
93  DAResidual& daResidual = const_cast<DAResidual&>(daResidual_);
94 
95  // zero all the matrices
96  MatZeroEntries(jacMat);
97 
98  // initialize residual vectors
99  Vec resVecRef, resVec;
100  VecDuplicate(wVec, &resVec);
101  VecDuplicate(wVec, &resVecRef);
102  VecZeroEntries(resVec);
103  VecZeroEntries(resVecRef);
104 
105  dictionary mOptions;
106  mOptions.set("updateState", 1);
107  mOptions.set("updateMesh", 0);
108  mOptions.set("setResVec", 1);
109  mOptions.set("isPC", options.getLabel("isPC"));
110  daResidual.masterFunction(mOptions, xvVec, wVec, resVecRef);
111 
112  scalar delta = daOption_.getSubDictOption<scalar>("adjPartDerivFDStep", "ACTD");
113  scalar rDelta = 1.0 / delta;
114  PetscScalar rDeltaValue = 0.0;
115  assignValueCheckAD(rDeltaValue, rDelta);
116 
117  Vec xvVecNew;
118  VecDuplicate(xvVec, &xvVecNew);
119  VecZeroEntries(xvVecNew);
120 
121  dictionary& diskModelSubDict = const_cast<dictionary&>(
122  daOption_.getAllOptions().subDict("fvSource").subDict(actuatorName));
123 
124  scalarList center;
125  diskModelSubDict.readEntry<scalarList>("center", center);
126  scalar innerRadius = diskModelSubDict.getScalar("innerRadius");
127  scalar outerRadius = diskModelSubDict.getScalar("outerRadius");
128  scalar POD = diskModelSubDict.getScalar("POD");
129  scalar scale = diskModelSubDict.getScalar("scale");
130  scalar expM = diskModelSubDict.getScalar("expM");
131  scalar expN = diskModelSubDict.getScalar("expN");
132  scalar targetThrust = diskModelSubDict.getScalar("targetThrust");
133 
134  DAFvSource& fvSource = const_cast<DAFvSource&>(
135  mesh_.thisDb().lookupObject<DAFvSource>("DAFvSource"));
136 
137  for (label i = 0; i < nActDVs_; i++)
138  {
139 
140  // perturb ACTD
141  if (i == 0)
142  {
143  // perturb x
144  center[0] += delta;
145  diskModelSubDict.set("center", center);
146  }
147  else if (i == 1)
148  {
149  // perturb y
150  center[1] += delta;
151  diskModelSubDict.set("center", center);
152  }
153  else if (i == 2)
154  {
155  // perturb z
156  center[2] += delta;
157  diskModelSubDict.set("center", center);
158  }
159  else if (i == 3)
160  {
161  // perturb innerRadius
162  innerRadius += delta;
163  diskModelSubDict.set("innerRadius", innerRadius);
164  }
165  else if (i == 4)
166  {
167  // perturb outerRadius
168  outerRadius += delta;
169  diskModelSubDict.set("outerRadius", outerRadius);
170  }
171  else if (i == 5)
172  {
173  // perturb scale
174  scale += delta;
175  diskModelSubDict.set("scale", scale);
176  }
177  else if (i == 6)
178  {
179  // perturb POD
180  POD += delta;
181  diskModelSubDict.set("POD", POD);
182  }
183  else if (i == 7)
184  {
185  // perturb expM
186  expM += delta;
187  diskModelSubDict.set("expM", expM);
188  }
189  else if (i == 8)
190  {
191  // perturb expN
192  expN += delta;
193  diskModelSubDict.set("expN", expN);
194  }
195  else if (i == 9)
196  {
197  targetThrust += delta;
198  diskModelSubDict.set("targetThrust", targetThrust);
199  }
200 
201  // we need to synchronize the DAOption to actuatorDVs
202  fvSource.syncDAOptionToActuatorDVs();
203 
204  // Info<<daOption_.getAllOptions().subDict("fvSource")<<endl;
205 
206  // compute residual
207  daResidual.masterFunction(mOptions, xvVecNew, wVec, resVec);
208 
209  // reset perturbation
210  if (i == 0)
211  {
212  // reset x
213  center[0] -= delta;
214  diskModelSubDict.set("center", center);
215  }
216  else if (i == 1)
217  {
218  // reset y
219  center[1] -= delta;
220  diskModelSubDict.set("center", center);
221  }
222  else if (i == 2)
223  {
224  // reset z
225  center[2] -= delta;
226  diskModelSubDict.set("center", center);
227  }
228  else if (i == 3)
229  {
230  // reset innerRadius
231  innerRadius -= delta;
232  diskModelSubDict.set("innerRadius", innerRadius);
233  }
234  else if (i == 4)
235  {
236  // reset outerRadius
237  outerRadius -= delta;
238  diskModelSubDict.set("outerRadius", outerRadius);
239  }
240  else if (i == 5)
241  {
242  // reset scale
243  scale -= delta;
244  diskModelSubDict.set("scale", scale);
245  }
246  else if (i == 6)
247  {
248  // reset POD
249  POD -= delta;
250  diskModelSubDict.set("POD", POD);
251  }
252  else if (i == 7)
253  {
254  // reset expM
255  expM -= delta;
256  diskModelSubDict.set("expM", expM);
257  }
258  else if (i == 8)
259  {
260  // reset expN
261  expN -= delta;
262  diskModelSubDict.set("expN", expN);
263  }
264  else if (i == 9)
265  {
266  targetThrust -= delta;
267  diskModelSubDict.set("targetThrust", targetThrust);
268  }
269 
270  // we need to synchronize the DAOption to actuatorDVs
271  fvSource.syncDAOptionToActuatorDVs();
272 
273  // Info<<daOption_.getAllOptions().subDict("fvSource")<<endl;
274 
275  // compute residual partial using finite-difference
276  VecAXPY(resVec, -1.0, resVecRef);
277  VecScale(resVec, rDeltaValue);
278 
279  // assign resVec to jacMat
280  PetscInt Istart, Iend;
281  VecGetOwnershipRange(resVec, &Istart, &Iend);
282 
283  const PetscScalar* resVecArray;
284  VecGetArrayRead(resVec, &resVecArray);
285  for (label j = Istart; j < Iend; j++)
286  {
287  label relIdx = j - Istart;
288  PetscScalar val = resVecArray[relIdx];
289  MatSetValue(jacMat, j, i, val, INSERT_VALUES);
290  }
291  VecRestoreArrayRead(resVec, &resVecArray);
292  }
293 
294  label eTime = mesh_.time().elapsedClockTime();
295  Info << modelType_ << " ExecutionTime: " << eTime << " s" << endl;
296 
297  // call the master function again to reset the xvVec and wVec to OpenFOAM fields and points
298  daResidual.masterFunction(mOptions, xvVec, wVec, resVecRef);
299 
300  MatAssemblyBegin(jacMat, MAT_FINAL_ASSEMBLY);
301  MatAssemblyEnd(jacMat, MAT_FINAL_ASSEMBLY);
302 }
303 
304 } // End namespace Foam
305 
306 // ************************************************************************* //
Foam::DAPartDerivdRdACTD::DAPartDerivdRdACTD
DAPartDerivdRdACTD(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAJacCon &daJacCon, const DAResidual &daResidual)
Definition: DAPartDerivdRdACTD.C:19
Foam::DAPartDerivdRdACTD::initializePartDerivMat
virtual void initializePartDerivMat(const dictionary &options, Mat jacMat)
initialize partial derivative matrix
Definition: DAPartDerivdRdACTD.C:40
Foam::DAFvSource
Definition: DAFvSource.H:34
Foam::DAJacCon
Definition: DAJacCon.H:34
Foam::DAPartDerivdRdACTD::calcPartDerivMat
virtual void calcPartDerivMat(const dictionary &options, const Vec xvVec, const Vec wVec, Mat jacMat)
compute the partial derivative matrix
Definition: DAPartDerivdRdACTD.C:69
Foam::DAOption
Definition: DAOption.H:29
Foam::DAPartDeriv::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAPartDeriv.H:60
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAPartDerivdRdACTD::nActDVs_
label nActDVs_
number of design variables for actuator disk
Definition: DAPartDerivdRdACTD.H:33
Foam::DAPartDeriv::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAPartDeriv.H:54
fvSource
volScalarField & fvSource
Definition: createRefsHeatTransfer.H:7
Foam::DAPartDeriv::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAPartDeriv.H:51
Foam::DAOption::getAllOptions
const dictionary & getAllOptions() const
return a reference of allOptions_ dictionary
Definition: DAOption.H:56
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
DAPartDerivdRdACTD.H
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::DAResidual
Definition: DAResidual.H:35
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::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DAPartDeriv::modelType_
const word modelType_
the name of the jacCon matrix
Definition: DAPartDeriv.H:48
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::DAIndex::nLocalAdjointStates
label nLocalAdjointStates
number of local adjoint states (including all cells and faces)
Definition: DAIndex.H:145