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