DAObjFuncPower.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAObjFuncPower.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAObjFuncPower, 0);
16 addToRunTimeSelectionTable(DAObjFunc, DAObjFuncPower, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  const fvMesh& mesh,
21  const DAOption& daOption,
22  const DAModel& daModel,
23  const DAIndex& daIndex,
24  const DAResidual& daResidual,
25  const word objFuncName,
26  const word objFuncPart,
27  const dictionary& objFuncDict)
28  : DAObjFunc(
29  mesh,
30  daOption,
31  daModel,
32  daIndex,
33  daResidual,
34  objFuncName,
35  objFuncPart,
36  objFuncDict),
37  daTurb_(daModel.getDATurbulenceModel())
38 {
39 
40  // for computing moment, first read in some parameters from objFuncDict_
41  // these parameters are only for moment objective
42 
43  // Assign type, this is common for all objectives
44  objFuncDict_.readEntry<word>("type", objFuncType_);
45 
46  scalarList dir;
47  objFuncDict_.readEntry<scalarList>("axis", dir);
48  momentDir_[0] = dir[0];
49  momentDir_[1] = dir[1];
50  momentDir_[2] = dir[2];
51 
52  if (fabs(mag(momentDir_) - 1.0) > 1.0e-4)
53  {
54  FatalErrorIn(" ") << "the magnitude of the axis parameter in "
55  << objFuncName << " " << objFuncPart << " is not 1.0!"
56  << abort(FatalError);
57  }
58 
59  scalarList center;
60  objFuncDict_.readEntry<scalarList>("center", center);
61  momentCenter_[0] = center[0];
62  momentCenter_[1] = center[1];
63  momentCenter_[2] = center[2];
64 
65  objFuncDict_.readEntry<scalar>("scale", scale_);
66 
67  // setup the connectivity for moment, this is needed in Foam::DAJacCondFdW
68  // need to determine the name of pressure because some buoyant OF solver use
69  // p_rgh as pressure
70  word pName = "p";
71  if (mesh_.thisDb().foundObject<volScalarField>("p_rgh"))
72  {
73  pName = "p_rgh";
74  }
75 #ifdef IncompressibleFlow
76  // For incompressible flow, it depends on zero level
77  // of U, nut, and p, and one level of U
78  objFuncConInfo_ = {
79  {"U", "nut", pName}, // level 0
80  {"U"}}; // level 1
81 #endif
82 
83 #ifdef CompressibleFlow
84  // For compressible flow, it depends on zero level
85  // of U, nut, T, and p, and one level of U
86  objFuncConInfo_ = {
87  {"U", "nut", "T", pName}, // level 0
88  {"U"}}; // level 1
89 #endif
90 
91  // now replace nut with the corrected name for the selected turbulence model
92  daModel.correctModelStates(objFuncConInfo_[0]);
93 }
94 
97  const labelList& objFuncFaceSources,
98  const labelList& objFuncCellSources,
99  scalarList& objFuncFaceValues,
100  scalarList& objFuncCellValues,
101  scalar& objFuncValue)
102 {
103  /*
104  Description:
105  Calculate the moment which consist of pressure and viscous components
106  of force cross-producting the r vector wrt to a ref point
107  This code for computiong force is modified based on:
108  src/functionObjects/forcces/forces.C
109 
110  Input:
111  objFuncFaceSources: List of face source (index) for this objective
112 
113  objFuncCellSources: List of cell source (index) for this objective
114 
115  Output:
116  objFuncFaceValues: the discrete value of objective for each face source (index).
117  This will be used for computing df/dw in the adjoint.
118 
119  objFuncCellValues: the discrete value of objective on each cell source (index).
120  This will be used for computing df/dw in the adjoint.
121 
122  objFuncValue: the sum of objective, reduced across all processsors and scaled by "scale"
123  */
124 
125  // reload the scale, which may be needed for multipoint optimization
126  objFuncDict_.readEntry<scalar>("scale", scale_);
127 
128  // initialize faceValues to zero
129  forAll(objFuncFaceValues, idxI)
130  {
131  objFuncFaceValues[idxI] = 0.0;
132  }
133  // initialize objFunValue
134  objFuncValue = 0.0;
135 
136  const objectRegistry& db = mesh_.thisDb();
137  const volScalarField& p = db.lookupObject<volScalarField>("p");
138 
139  const surfaceVectorField::Boundary& Sfb = mesh_.Sf().boundaryField();
140 
141  tmp<volSymmTensorField> tdevRhoReff = daTurb_.devRhoReff();
142  const volSymmTensorField::Boundary& devRhoReffb = tdevRhoReff().boundaryField();
143 
144  // calculate discrete force for each objFuncFace
145  forAll(objFuncFaceSources, idxI)
146  {
147  const label& objFuncFaceI = objFuncFaceSources[idxI];
148  label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
149  const label patchI = daIndex_.bFacePatchI[bFaceI];
150  const label faceI = daIndex_.bFaceFaceI[bFaceI];
151 
152  // normal force
153  vector fN(Sfb[patchI][faceI] * p.boundaryField()[patchI][faceI]);
154  // tangential force
155  vector fT(Sfb[patchI][faceI] & devRhoReffb[patchI][faceI]);
156  // r vector
157  vector rVec = mesh_.Cf().boundaryField()[patchI][faceI] - momentCenter_;
158  // compute moment
159  objFuncFaceValues[idxI] = scale_ * (rVec ^ (fN + fT)) & momentDir_;
160 
161  objFuncValue += objFuncFaceValues[idxI];
162  }
163 
164  const IOMRFZoneListDF& MRF = mesh_.thisDb().lookupObject<IOMRFZoneListDF>("MRFProperties");
165  scalar& omega = const_cast<scalar&>(MRF.getOmegaRef());
166  objFuncValue *= omega;
167 
168  // need to reduce the sum of force across all processors
169  reduce(objFuncValue, sumOp<scalar>());
170 
171  return;
172 }
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 } // End namespace Foam
176 
177 // ************************************************************************* //
Foam::DAObjFunc::objFuncType_
word objFuncType_
the type of the objective function
Definition: DAObjFunc.H:66
forAll
forAll(pseudoP.boundaryField(), patchI)
Definition: solvePseudoPEqn.H:10
Foam::DAIndex::bFacePatchI
labelList bFacePatchI
Definition: DAIndex.H:110
Foam::DAOption
Definition: DAOption.H:29
Foam::DAObjFuncPower::momentCenter_
vector momentCenter_
the center of rotation for moment
Definition: DAObjFuncPower.H:36
Foam::DAObjFunc::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAObjFunc.H:54
daOption
DAOption daOption(mesh, pyOptions_)
MRF
IOMRFZoneListDF & MRF
Definition: createRefsRhoSimple.H:18
p
volScalarField & p
Definition: createRefsPimple.H:6
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAObjFunc::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAObjFunc.H:45
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAModel
Definition: DAModel.H:59
Foam::DAObjFuncPower::DAObjFuncPower
DAObjFuncPower(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: DAObjFuncPower.C:19
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAObjFuncPower::momentDir_
vector momentDir_
the direction of the moment
Definition: DAObjFuncPower.H:33
Foam::DATurbulenceModel::devRhoReff
tmp< volSymmTensorField > devRhoReff() const
dev terms
Definition: DATurbulenceModel.C:306
Foam::IOMRFZoneListDF
Definition: IOMRFZoneListDF.H:49
Foam::DAObjFuncPower::calcObjFunc
virtual void calcObjFunc(const labelList &objFuncFaceSources, const labelList &objFuncCellSources, scalarList &objFuncFaceValues, scalarList &objFuncCellValues, scalar &objFuncValue)
calculate the value of objective function
Definition: DAObjFuncPower.C:96
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DAResidual
Definition: DAResidual.H:35
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DAObjFunc::scale_
scalar scale_
scale of the objective function
Definition: DAObjFunc.H:87
daModel
DAModel daModel(mesh, daOption)
Foam::DAObjFunc::objFuncConInfo_
List< List< word > > objFuncConInfo_
the connectivity information for the objective function
Definition: DAObjFunc.H:93
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAObjFuncPower::daTurb_
const DATurbulenceModel & daTurb_
DATurbulenceModel object.
Definition: DAObjFuncPower.H:39
Foam::DAObjFunc::objFuncDict_
const dictionary & objFuncDict_
dictionary containing the information for the objective function
Definition: DAObjFunc.H:69
Foam::DAIndex::bFaceFaceI
labelList bFaceFaceI
Definition: DAIndex.H:116
Foam::DAObjFunc
Definition: DAObjFunc.H:33
DAObjFuncPower.H