DAObjFuncPatchMean.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAObjFuncPatchMean.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAObjFuncPatchMean, 0);
16 addToRunTimeSelectionTable(DAObjFunc, DAObjFuncPatchMean, 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 {
38  // Assign type, this is common for all objectives
39  objFuncDict_.readEntry<word>("type", objFuncType_);
40 
41  objFuncDict_.readEntry<scalar>("scale", scale_);
42 
43  objFuncDict_.readEntry<word>("varName", varName_);
44 
45  objFuncDict_.readEntry<word>("varType", varType_);
46 
47  objFuncDict_.readEntry<label>("component", component_);
48 }
49 
52  const labelList& objFuncFaceSources,
53  const labelList& objFuncCellSources,
54  scalarList& objFuncFaceValues,
55  scalarList& objFuncCellValues,
56  scalar& objFuncValue)
57 {
58  /*
59  Description:
60  Calculate the patch mean
61 
62  Input:
63  objFuncFaceSources: List of face source (index) for this objective
64 
65  objFuncCellSources: List of cell source (index) for this objective
66 
67  Output:
68  objFuncFaceValues: the discrete value of objective for each face source (index).
69  This will be used for computing df/dw in the adjoint.
70 
71  objFuncCellValues: the discrete value of objective on each cell source (index).
72  This will be used for computing df/dw in the adjoint.
73 
74  objFuncValue: the sum of objective, reduced across all processors and scaled by "scale"
75  */
76 
77  // calculate the area of all the patches. We need to recompute because the surface area
78  // may change during the optimization
79  areaSum_ = 0.0;
80  forAll(objFuncFaceSources, idxI)
81  {
82  const label& objFuncFaceI = objFuncFaceSources[idxI];
83  label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
84  const label patchI = daIndex_.bFacePatchI[bFaceI];
85  const label faceI = daIndex_.bFaceFaceI[bFaceI];
86  areaSum_ += mesh_.magSf().boundaryField()[patchI][faceI];
87  }
88  reduce(areaSum_, sumOp<scalar>());
89 
90  // initialize objFunValue
91  objFuncValue = 0.0;
92 
93  const objectRegistry& db = mesh_.thisDb();
94 
95  if (varType_ == "scalar")
96  {
97  const volScalarField& var = db.lookupObject<volScalarField>(varName_);
98 
99  // calculate area weighted heat flux
100  forAll(objFuncFaceSources, idxI)
101  {
102  const label& objFuncFaceI = objFuncFaceSources[idxI];
103  label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
104  const label patchI = daIndex_.bFacePatchI[bFaceI];
105  const label faceI = daIndex_.bFaceFaceI[bFaceI];
106  scalar area = mesh_.magSf().boundaryField()[patchI][faceI];
107  objFuncValue += scale_ * area * var.boundaryField()[patchI][faceI] / areaSum_;
108  }
109  }
110  else if (varType_ == "vector")
111  {
112  const volVectorField& var = db.lookupObject<volVectorField>(varName_);
113 
114  // calculate area weighted heat flux
115  forAll(objFuncFaceSources, idxI)
116  {
117  const label& objFuncFaceI = objFuncFaceSources[idxI];
118  label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
119  const label patchI = daIndex_.bFacePatchI[bFaceI];
120  const label faceI = daIndex_.bFaceFaceI[bFaceI];
121  scalar area = mesh_.magSf().boundaryField()[patchI][faceI];
122  objFuncValue += scale_ * area * var.boundaryField()[patchI][faceI][component_] / areaSum_;
123  }
124  }
125  else
126  {
127  FatalErrorIn("DAObjFuncPatchMean::calcObjFunc")
128  << "varType not valid. Options are scalar or vector"
129  << abort(FatalError);
130  }
131 
132  // need to reduce the sum of force across all processors
133  reduce(objFuncValue, sumOp<scalar>());
134 
135  return;
136 }
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 } // End namespace Foam
141 
142 // ************************************************************************* //
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::DAObjFuncPatchMean::component_
label component_
if vector which element?
Definition: DAObjFuncPatchMean.H:42
Foam::DAObjFuncPatchMean::calcObjFunc
virtual void calcObjFunc(const labelList &objFuncFaceSources, const labelList &objFuncCellSources, scalarList &objFuncFaceValues, scalarList &objFuncCellValues, scalar &objFuncValue)
calculate the value of objective function
Definition: DAObjFuncPatchMean.C:51
Foam::DAIndex::bFacePatchI
labelList bFacePatchI
Definition: DAIndex.H:110
Foam::DAOption
Definition: DAOption.H:29
Foam::DAObjFunc::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAObjFunc.H:54
daOption
DAOption daOption(mesh, pyOptions_)
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::DAObjFuncPatchMean::varName_
word varName_
name of the variable
Definition: DAObjFuncPatchMean.H:36
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DAObjFuncPatchMean::areaSum_
scalar areaSum_
the area of all total pressure patches
Definition: DAObjFuncPatchMean.H:33
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::DAObjFuncPatchMean::varType_
word varType_
type of the variable either vector or scalar
Definition: DAObjFuncPatchMean.H:39
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAObjFunc::objFuncDict_
const dictionary & objFuncDict_
dictionary containing the information for the objective function
Definition: DAObjFunc.H:69
DAObjFuncPatchMean.H
Foam::DAIndex::bFaceFaceI
labelList bFaceFaceI
Definition: DAIndex.H:116
Foam::DAObjFuncPatchMean::DAObjFuncPatchMean
DAObjFuncPatchMean(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: DAObjFuncPatchMean.C:19
Foam::DAObjFunc
Definition: DAObjFunc.H:33