DAPartDerivdRdACTL.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAPartDerivdRdACTL.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAPartDerivdRdACTL, 0);
16 addToRunTimeSelectionTable(DAPartDeriv, DAPartDerivdRdACTL, 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 dRdACTLT
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 dRdACTL 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", "ACTL");
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& pointModelSubDict = const_cast<dictionary&>(
122  daOption_.getAllOptions().subDict("fvSource").subDict(actuatorName));
123 
124  scalarList center;
125  pointModelSubDict.readEntry<scalarList>("center", center);
126  scalar innerRadius = pointModelSubDict.getScalar("innerRadius");
127  scalar outerRadius = pointModelSubDict.getScalar("outerRadius");
128  scalar rpm = pointModelSubDict.getScalar("rpm");
129  scalar phase = pointModelSubDict.getScalar("phase");
130  scalar scale = pointModelSubDict.getScalar("scale");
131  scalar POD = pointModelSubDict.getScalar("POD");
132  scalar expM = pointModelSubDict.getScalar("expM");
133  scalar expN = pointModelSubDict.getScalar("expN");
134 
135  for (label i = 0; i < nActDVs_; i++)
136  {
137 
138  // perturb ACTL
139  if (i == 0)
140  {
141  // perturb x
142  center[0] += delta;
143  pointModelSubDict.set("center", center);
144  }
145  else if (i == 1)
146  {
147  // perturb y
148  center[1] += delta;
149  pointModelSubDict.set("center", center);
150  }
151  else if (i == 2)
152  {
153  // perturb z
154  center[2] += delta;
155  pointModelSubDict.set("center", center);
156  }
157  else if (i == 3)
158  {
159  // perturb innerRadius
160  innerRadius += delta;
161  pointModelSubDict.set("innerRadius", innerRadius);
162  }
163  else if (i == 4)
164  {
165  // perturb outerRadius
166  outerRadius += delta;
167  pointModelSubDict.set("outerRadius", outerRadius);
168  }
169  else if (i == 5)
170  {
171  // perturb rpm
172  rpm += delta;
173  pointModelSubDict.set("rpm", rpm);
174  }
175  else if (i == 6)
176  {
177  // perturb phase
178  phase += delta;
179  pointModelSubDict.set("phase", phase);
180  }
181  else if (i == 7)
182  {
183  // perturb scale
184  scale += delta;
185  pointModelSubDict.set("scale", scale);
186  }
187  else if (i == 8)
188  {
189  // perturb POD
190  POD += delta;
191  pointModelSubDict.set("POD", POD);
192  }
193  else if (i == 9)
194  {
195  // perturb expM
196  expM += delta;
197  pointModelSubDict.set("expM", expM);
198  }
199  else if (i == 10)
200  {
201  // perturb expN
202  expN += delta;
203  pointModelSubDict.set("expN", expN);
204  }
205 
206  // Info<<daOption_.getAllOptions().subDict("fvSource")<<endl;
207 
208  // compute residual
209  daResidual.masterFunction(mOptions, xvVecNew, wVec, resVec);
210 
211  // reset perturbation
212  if (i == 0)
213  {
214  // reset x
215  center[0] -= delta;
216  pointModelSubDict.set("center", center);
217  }
218  else if (i == 1)
219  {
220  // reset y
221  center[1] -= delta;
222  pointModelSubDict.set("center", center);
223  }
224  else if (i == 2)
225  {
226  // reset z
227  center[2] -= delta;
228  pointModelSubDict.set("center", center);
229  }
230  else if (i == 3)
231  {
232  // reset innerRadius
233  innerRadius -= delta;
234  pointModelSubDict.set("innerRadius", innerRadius);
235  }
236  else if (i == 4)
237  {
238  // reset outerRadius
239  outerRadius -= delta;
240  pointModelSubDict.set("outerRadius", outerRadius);
241  }
242  else if (i == 5)
243  {
244  // reset rpm
245  rpm -= delta;
246  pointModelSubDict.set("rpm", rpm);
247  }
248  else if (i == 6)
249  {
250  // reset phase
251  phase -= delta;
252  pointModelSubDict.set("phase", phase);
253  }
254  else if (i == 7)
255  {
256  // reset scale
257  scale -= delta;
258  pointModelSubDict.set("scale", scale);
259  }
260  else if (i == 8)
261  {
262  // reset POD
263  POD -= delta;
264  pointModelSubDict.set("POD", POD);
265  }
266  else if (i == 9)
267  {
268  // reset expM
269  expM -= delta;
270  pointModelSubDict.set("expM", expM);
271  }
272  else if (i == 10)
273  {
274  // reset expN
275  expN -= delta;
276  pointModelSubDict.set("expN", expN);
277  }
278 
279  // Info<<daOption_.getAllOptions().subDict("fvSource")<<endl;
280 
281  // compute residual partial using finite-difference
282  VecAXPY(resVec, -1.0, resVecRef);
283  VecScale(resVec, rDeltaValue);
284 
285  // assign resVec to jacMat
286  PetscInt Istart, Iend;
287  VecGetOwnershipRange(resVec, &Istart, &Iend);
288 
289  const PetscScalar* resVecArray;
290  VecGetArrayRead(resVec, &resVecArray);
291  for (label j = Istart; j < Iend; j++)
292  {
293  label relIdx = j - Istart;
294  PetscScalar val = resVecArray[relIdx];
295  MatSetValue(jacMat, j, i, val, INSERT_VALUES);
296  }
297  VecRestoreArrayRead(resVec, &resVecArray);
298  }
299 
300  label eTime = mesh_.time().elapsedClockTime();
301  Info << modelType_ << " ExecutionTime: " << eTime << " s" << endl;
302 
303  // call the master function again to reset the xvVec and wVec to OpenFOAM fields and points
304  daResidual.masterFunction(mOptions, xvVec, wVec, resVecRef);
305 
306  MatAssemblyBegin(jacMat, MAT_FINAL_ASSEMBLY);
307  MatAssemblyEnd(jacMat, MAT_FINAL_ASSEMBLY);
308 }
309 
310 } // End namespace Foam
311 
312 // ************************************************************************* //
Foam::DAPartDerivdRdACTL::calcPartDerivMat
virtual void calcPartDerivMat(const dictionary &options, const Vec xvVec, const Vec wVec, Mat jacMat)
compute the partial derivative matrix
Definition: DAPartDerivdRdACTL.C:69
Foam::DAJacCon
Definition: DAJacCon.H:34
Foam::DAPartDerivdRdACTL::initializePartDerivMat
virtual void initializePartDerivMat(const dictionary &options, Mat jacMat)
initialize partial derivative matrix
Definition: DAPartDerivdRdACTL.C:40
Foam::DAOption
Definition: DAOption.H:29
Foam::DAPartDeriv::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAPartDeriv.H:60
daOption
DAOption daOption(mesh, pyOptions_)
DAPartDerivdRdACTL.H
Foam::DAPartDeriv::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAPartDeriv.H:54
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
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::DAPartDerivdRdACTL::nActDVs_
label nActDVs_
number of design variables for actuator points
Definition: DAPartDerivdRdACTL.H:33
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::DAPartDerivdRdACTL::DAPartDerivdRdACTL
DAPartDerivdRdACTL(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAJacCon &daJacCon, const DAResidual &daResidual)
Definition: DAPartDerivdRdACTL.C:19
Foam::DAIndex::nLocalAdjointStates
label nLocalAdjointStates
number of local adjoint states (including all cells and faces)
Definition: DAIndex.H:145