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