DAPartDerivdRdW.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAPartDerivdRdW.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAPartDerivdRdW, 0);
16 addToRunTimeSelectionTable(DAPartDeriv, DAPartDerivdRdW, 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.transposed. Whether to compute the transposed of dRdW
50  */
51 
52  label transposed = options.getLabel("transposed");
53 
54  // now initialize the memory for the jacobian itself
55  label localSize = daIndex_.nLocalAdjointStates;
56 
57  // create dRdWT
58  //MatCreate(PETSC_COMM_WORLD, jacMat);
59  MatSetSizes(
60  jacMat,
61  localSize,
62  localSize,
63  PETSC_DETERMINE,
64  PETSC_DETERMINE);
65  MatSetFromOptions(jacMat);
66  daJacCon_.preallocatedRdW(jacMat, transposed);
67  //MatSetOption(jacMat, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
68  MatSetUp(jacMat);
69  MatZeroEntries(jacMat);
70  Info << "Partial derivative matrix created. " << mesh_.time().elapsedClockTime() << " s" << endl;
71 }
72 
74  const dictionary& options,
75  const Vec xvVec,
76  const Vec wVec,
77  Mat jacMat)
78 {
79  /*
80  Description:
81  Compute jacMat. We use coloring accelerated finite-difference
82 
83  Input:
84 
85  options.transposed. Whether to compute the transposed of dRdW
86 
87  options.isPC: whether to compute the jacMat for preconditioner
88 
89  options.lowerBound: any |value| that is smaller than lowerBound will be set to zero in dRdW
90 
91  xvVec: the volume mesh coordinate vector
92 
93  wVec: the state variable vector
94 
95  Output:
96  jacMat: the partial derivative matrix dRdW to compute
97  */
98 
99  label transposed = options.getLabel("transposed");
100 
101  // initialize coloredColumn vector
102  Vec coloredColumn;
103  VecDuplicate(wVec, &coloredColumn);
104  VecZeroEntries(coloredColumn);
105 
106  DAResidual& daResidual = const_cast<DAResidual&>(daResidual_);
107 
108  // zero all the matrices
109  MatZeroEntries(jacMat);
110 
111  Vec wVecNew;
112  VecDuplicate(wVec, &wVecNew);
113  VecCopy(wVec, wVecNew);
114 
115  // initialize residual vectors
116  Vec resVecRef, resVec;
117  VecDuplicate(wVec, &resVec);
118  VecDuplicate(wVec, &resVecRef);
119  VecZeroEntries(resVec);
120  VecZeroEntries(resVecRef);
121 
122  // set up state normalization vector
123  Vec normStatePerturbVec;
124  this->setNormStatePerturbVec(&normStatePerturbVec);
125 
126  dictionary mOptions;
127  mOptions.set("updateState", 1);
128  mOptions.set("updateMesh", 0);
129  mOptions.set("setResVec", 1);
130  mOptions.set("isPC", options.getLabel("isPC"));
131  daResidual.masterFunction(mOptions, xvVec, wVec, resVecRef);
132 
133  scalar jacLowerBound = options.getScalar("lowerBound");
134 
135  scalar delta = daOption_.getSubDictOption<scalar>("adjPartDerivFDStep", "State");
136  scalar rDelta = 1.0 / delta;
137  PetscScalar rDeltaValue = 0.0;
138  assignValueCheckAD(rDeltaValue, rDelta);
139 
140  label nColors = daJacCon_.getNJacConColors();
141 
142  word partDerivName = modelType_;
143  if (transposed)
144  {
145  partDerivName += "T";
146  }
147  if (options.getLabel("isPC"))
148  {
149  partDerivName += "PC";
150  }
151 
152  label printInterval = daOption_.getOption<label>("printInterval");
153  for (label color = 0; color < nColors; color++)
154  {
155  label eTime = mesh_.time().elapsedClockTime();
156  // print progress
157  if (color % printInterval == 0 or color == nColors - 1)
158  {
159  Info << partDerivName << ": " << color << " of " << nColors
160  << ", ExecutionTime: " << eTime << " s" << endl;
161  }
162 
163  // perturb states
164  this->perturbStates(
166  normStatePerturbVec,
167  color,
168  delta,
169  wVecNew);
170 
171  // compute residual
172  daResidual.masterFunction(mOptions, xvVec, wVecNew, resVec);
173 
174  // reset state perburbation
175  VecCopy(wVec, wVecNew);
176 
177  // compute residual partial using finite-difference
178  VecAXPY(resVec, -1.0, resVecRef);
179  VecScale(resVec, rDeltaValue);
180 
181  // compute the colored coloumn and assign resVec to jacMat
182  daJacCon_.calcColoredColumns(color, coloredColumn);
183  this->setPartDerivMat(resVec, coloredColumn, transposed, jacMat, jacLowerBound);
184  }
185 
186  // call masterFunction again to reset the wVec to OpenFOAM field
187  daResidual.masterFunction(mOptions, xvVec, wVec, resVecRef);
188 
189  MatAssemblyBegin(jacMat, MAT_FINAL_ASSEMBLY);
190  MatAssemblyEnd(jacMat, MAT_FINAL_ASSEMBLY);
191 
192  if (daOption_.getOption<label>("debug"))
193  {
194  daIndex_.printMatChars(jacMat);
195  }
196 }
197 
198 } // End namespace Foam
199 
200 // ************************************************************************* //
Foam::DAJacCon
Definition: DAJacCon.H:34
Foam::DAPartDerivdRdW::calcPartDerivMat
virtual void calcPartDerivMat(const dictionary &options, const Vec xvVec, const Vec wVec, Mat jacMat)
compute the partial derivative matrix
Definition: DAPartDerivdRdW.C:73
Foam::DAJacCon::calcColoredColumns
void calcColoredColumns(const label colorI, Vec coloredColumn) const
calculate the colored column vector
Definition: DAJacCon.C:1950
Foam::DAOption
Definition: DAOption.H:29
Foam::DAPartDeriv::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAPartDeriv.H:60
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAPartDeriv::setPartDerivMat
void setPartDerivMat(const Vec resVec, const Vec coloredColumn, const label transposed, Mat jacMat, const scalar jacLowerBound=1e-30) const
set values for the partial derivative matrix
Definition: DAPartDeriv.C:165
Foam::DAPartDeriv::daJacCon_
const DAJacCon & daJacCon_
DAJacCon object.
Definition: DAPartDeriv.H:63
Foam::DAJacCon::getJacConColor
Vec getJacConColor() const
return DAJacCon::jacConColors_
Definition: DAJacCon.H:272
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::DAJacCon::getNJacConColors
label getNJacConColors() const
get the number of JacCon colors
Definition: DAJacCon.H:252
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::DAPartDerivdRdW::initializePartDerivMat
virtual void initializePartDerivMat(const dictionary &options, Mat jacMat)
initialize partial derivative matrix
Definition: DAPartDerivdRdW.C:40
Foam::DAIndex
Definition: DAIndex.H:32
assignValueCheckAD
#define assignValueCheckAD(valA, valB)
Definition: DAMacroFunctions.H:103
Foam::DAJacCon::preallocatedRdW
virtual void preallocatedRdW(Mat dRMat, const label transposed) const
preallocate dRdW matrix using the preallocVec
Definition: DAJacCon.C:1915
Foam::DAModel
Definition: DAModel.H:59
DAPartDerivdRdW.H
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::DAIndex::printMatChars
void printMatChars(const Mat matIn) const
print petsc matrix statistics such as the max on/off diagonral ratio
Definition: DAIndex.C:1010
Foam::DAPartDeriv::daResidual_
const DAResidual & daResidual_
DAResidual object.
Definition: DAPartDeriv.H:66
Foam::DAPartDeriv::perturbStates
void perturbStates(const Vec jacConColors, const Vec normStatePerturbVec, const label colorI, const scalar delta, Vec wVec)
perturb state variables given a color index
Definition: DAPartDeriv.C:98
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAPartDeriv::setNormStatePerturbVec
void setNormStatePerturbVec(Vec *normStatePerturbVec)
setup the state normalization vector
Definition: DAPartDeriv.C:511
Foam::DAPartDerivdRdW::DAPartDerivdRdW
DAPartDerivdRdW(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const DAJacCon &daJacCon, const DAResidual &daResidual)
Definition: DAPartDerivdRdW.C:19
Foam::DAIndex::nLocalAdjointStates
label nLocalAdjointStates
number of local adjoint states (including all cells and faces)
Definition: DAIndex.H:145