DAPartDerivdRdFFD.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAPartDerivdRdFFD.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAPartDerivdRdFFD, 0);
16 addToRunTimeSelectionTable(DAPartDeriv, DAPartDerivdRdFFD, 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.nDesignVars. The number of design variable for dRdFFD
50  */
51 
52  label nDesignVars = options.getLabel("nDesignVars");
53 
54  // now initialize the memory for the jacobian itself
55  label localSize = daIndex_.nLocalAdjointStates;
56 
57  // create dRdFFDT
58  //MatCreate(PETSC_COMM_WORLD, jacMat);
59  MatSetSizes(
60  jacMat,
61  localSize,
62  PETSC_DECIDE,
63  PETSC_DETERMINE,
64  nDesignVars);
65  MatSetFromOptions(jacMat);
66  MatMPIAIJSetPreallocation(jacMat, nDesignVars, NULL, nDesignVars, NULL);
67  MatSeqAIJSetPreallocation(jacMat, nDesignVars, NULL);
68  //MatSetOption(jacMat, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
69  MatSetUp(jacMat);
70  MatZeroEntries(jacMat);
71  Info << "Partial derivative matrix created. " << mesh_.time().elapsedClockTime() << " s" << endl;
72 }
73 
75  const dictionary& options,
76  const Vec xvVec,
77  const Vec wVec,
78  Mat jacMat)
79 {
80  /*
81  Description:
82  Compute jacMat. We use brute-force finite-difference
83 
84  Input:
85 
86  options.nDesignVars. The number of design variable for dRdFFD
87 
88  options.isPC: whether to compute the jacMat for preconditioner
89 
90  xvVec: the volume mesh coordinate vector
91 
92  wVec: the state variable vector
93 
94  Output:
95  jacMat: the partial derivative matrix dRdFFD to compute
96  */
97 
98  label nDesignVars = options.getLabel("nDesignVars");
99 
100  DAResidual& daResidual = const_cast<DAResidual&>(daResidual_);
101 
102  // zero all the matrices
103  MatZeroEntries(jacMat);
104 
105  // initialize residual vectors
106  Vec resVecRef, resVec;
107  VecDuplicate(wVec, &resVec);
108  VecDuplicate(wVec, &resVecRef);
109  VecZeroEntries(resVec);
110  VecZeroEntries(resVecRef);
111 
112  dictionary mOptions;
113  mOptions.set("updateState", 1);
114  mOptions.set("updateMesh", 1);
115  mOptions.set("setResVec", 1);
116  mOptions.set("isPC", options.getLabel("isPC"));
117  daResidual.masterFunction(mOptions, xvVec, wVec, resVecRef);
118 
119  scalar delta = daOption_.getSubDictOption<scalar>("adjPartDerivFDStep", "FFD");
120  scalar rDelta = 1.0 / delta;
121  PetscScalar rDeltaValue = 0.0;
122  assignValueCheckAD(rDeltaValue, rDelta);
123 
124  Vec xvVecNew;
125  VecDuplicate(xvVec, &xvVecNew);
126  VecZeroEntries(xvVecNew);
127 
128  label printInterval = daOption_.getOption<label>("printInterval");
129  for (label i = 0; i < nDesignVars; i++)
130  {
131  label eTime = mesh_.time().elapsedClockTime();
132  // print progress
133  if (i % printInterval == 0 or i == nDesignVars - 1)
134  {
135  Info << modelType_ << ": " << i << " of " << nDesignVars
136  << ", ExecutionTime: " << eTime << " s" << endl;
137  }
138 
139  // perturb FFD
140  VecZeroEntries(xvVecNew);
141  MatGetColumnVector(dXvdFFDMat_, xvVecNew, i);
142  VecAXPY(xvVecNew, 1.0, xvVec);
143 
144  // compute residual
145  daResidual.masterFunction(mOptions, xvVecNew, wVec, resVec);
146 
147  // no need to reset FFD here
148 
149  // compute residual partial using finite-difference
150  VecAXPY(resVec, -1.0, resVecRef);
151  VecScale(resVec, rDeltaValue);
152 
153  // assign resVec to jacMat
154  PetscInt Istart, Iend;
155  VecGetOwnershipRange(resVec, &Istart, &Iend);
156 
157  const PetscScalar* resVecArray;
158  VecGetArrayRead(resVec, &resVecArray);
159  for (label j = Istart; j < Iend; j++)
160  {
161  label relIdx = j - Istart;
162  PetscScalar val = resVecArray[relIdx];
163  MatSetValue(jacMat, j, i, val, INSERT_VALUES);
164  }
165  VecRestoreArrayRead(resVec, &resVecArray);
166  }
167 
168  // call the master function again to reset the xvVec and wVec to OpenFOAM fields and points
169  daResidual.masterFunction(mOptions, xvVec, wVec, resVecRef);
170 
171  MatAssemblyBegin(jacMat, MAT_FINAL_ASSEMBLY);
172  MatAssemblyEnd(jacMat, MAT_FINAL_ASSEMBLY);
173 }
174 
175 } // End namespace Foam
176 
177 // ************************************************************************* //
Foam::DAPartDerivdRdFFD::initializePartDerivMat
virtual void initializePartDerivMat(const dictionary &options, Mat jacMat)
initialize partial derivative matrix
Definition: DAPartDerivdRdFFD.C:40
Foam::DAJacCon
Definition: DAJacCon.H:34
Foam::DAPartDerivdRdFFD::calcPartDerivMat
virtual void calcPartDerivMat(const dictionary &options, const Vec xvVec, const Vec wVec, Mat jacMat)
compute the partial derivative matrix
Definition: DAPartDerivdRdFFD.C:74
Foam::DAOption
Definition: DAOption.H:29
DAPartDerivdRdFFD.H
Foam::DAPartDerivdRdFFD::DAPartDerivdRdFFD
DAPartDerivdRdFFD(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAJacCon &daJacCon, const DAResidual &daResidual)
Definition: DAPartDerivdRdFFD.C:19
Foam::DAPartDeriv::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAPartDeriv.H:60
daOption
DAOption daOption(mesh, pyOptions_)
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::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::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
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::DAIndex::nLocalAdjointStates
label nLocalAdjointStates
number of local adjoint states (including all cells and faces)
Definition: DAIndex.H:145