DAObjFuncTotalPressure.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAObjFuncTotalPressure, 0);
16 addToRunTimeSelectionTable(DAObjFunc, DAObjFuncTotalPressure, 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_(const_cast<DATurbulenceModel&>(daModel.getDATurbulenceModel())),
38  rho_(daTurb_.getRho())
39 {
40  // Assign type, this is common for all objectives
41  objFuncDict_.readEntry<word>("type", objFuncType_);
42 
43  // setup the connectivity for total pressure, this is needed in Foam::DAJacCondFdW
44  // NOTE: for pressureInlet velocity, U depends on phi
45 #ifdef CompressibleFlow
46  // for compressible cases
47  objFuncConInfo_ = {
48  {"U", "T", "p", "phi"}, // level 0
49  {"U", "T", "p", "phi"}}; // level 1
50 #endif
51 #ifdef IncompressibleFlow
52  // for incompressible cases
53  objFuncConInfo_ = {
54  {"U", "p", "phi"}, // level 0
55  {"U", "p", "phi"}}; // level 1
56 #endif
57 
58  objFuncDict_.readEntry<scalar>("scale", scale_);
59 }
60 
63  const labelList& objFuncFaceSources,
64  const labelList& objFuncCellSources,
65  scalarList& objFuncFaceValues,
66  scalarList& objFuncCellValues,
67  scalar& objFuncValue)
68 {
69  /*
70  Description:
71  Calculate the total pressure TP=p+0.5*rho*U^2.
72 
73  Input:
74  objFuncFaceSources: List of face source (index) for this objective
75 
76  objFuncCellSources: List of cell source (index) for this objective
77 
78  Output:
79  objFuncFaceValues: the discrete value of objective for each face source (index).
80  This will be used for computing df/dw in the adjoint.
81 
82  objFuncCellValues: the discrete value of objective on each cell source (index).
83  This will be used for computing df/dw in the adjoint.
84 
85  objFuncValue: the sum of objective, reduced across all processors and scaled by "scale"
86  */
87 
88  // calculate the area of all the heat flux patches
89  if (areaSum_ < 0.0)
90  {
91  areaSum_ = 0.0;
92  forAll(objFuncFaceSources, idxI)
93  {
94  const label& objFuncFaceI = objFuncFaceSources[idxI];
95  label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
96  const label patchI = daIndex_.bFacePatchI[bFaceI];
97  const label faceI = daIndex_.bFaceFaceI[bFaceI];
98  areaSum_ += mesh_.magSf().boundaryField()[patchI][faceI];
99  }
100  reduce(areaSum_, sumOp<scalar>());
101  }
102 
103  // initialize faceValues to zero
104  forAll(objFuncFaceValues, idxI)
105  {
106  objFuncFaceValues[idxI] = 0.0;
107  }
108  // initialize objFunValue
109  objFuncValue = 0.0;
110 
111  const objectRegistry& db = mesh_.thisDb();
112  const volScalarField& p = db.lookupObject<volScalarField>("p");
113  const volVectorField& U = db.lookupObject<volVectorField>("U");
114 
115  const volScalarField::Boundary& pBf = p.boundaryField();
116  const volVectorField::Boundary& UBf = U.boundaryField();
117  const volScalarField::Boundary& rhoBf = rho_.boundaryField();
118 
119  // calculate area weighted heat flux
120  forAll(objFuncFaceSources, idxI)
121  {
122  const label& objFuncFaceI = objFuncFaceSources[idxI];
123  label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
124  const label patchI = daIndex_.bFacePatchI[bFaceI];
125  const label faceI = daIndex_.bFaceFaceI[bFaceI];
126 
127  scalar area = mesh_.magSf().boundaryField()[patchI][faceI];
128  objFuncFaceValues[idxI] =
129  pBf[patchI][faceI] + 0.5 * rhoBf[patchI][faceI] * mag(UBf[patchI][faceI]) * mag(UBf[patchI][faceI]);
130  objFuncFaceValues[idxI] *= scale_ * area / areaSum_;
131 
132  objFuncValue += objFuncFaceValues[idxI];
133  }
134 
135  // need to reduce the sum of force across all processors
136  reduce(objFuncValue, sumOp<scalar>());
137 
138  return;
139 }
140 
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
142 
143 } // End namespace Foam
144 
145 // ************************************************************************* //
Foam::DAObjFunc::objFuncType_
word objFuncType_
the type of the objective function
Definition: DAObjFunc.H:66
U
U
Definition: pEqnPimpleDyM.H:60
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::DAObjFunc::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAObjFunc.H:54
daOption
DAOption daOption(mesh, pyOptions_)
p
volScalarField & p
Definition: createRefsPimple.H:6
DAObjFuncTotalPressure.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAObjFunc::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAObjFunc.H:45
Foam::DAObjFuncTotalPressure::rho_
volScalarField rho_
density
Definition: DAObjFuncTotalPressure.H:36
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAObjFuncTotalPressure::calcObjFunc
virtual void calcObjFunc(const labelList &objFuncFaceSources, const labelList &objFuncCellSources, scalarList &objFuncFaceValues, scalarList &objFuncCellValues, scalar &objFuncValue)
calculate the value of objective function
Definition: DAObjFuncTotalPressure.C:62
Foam::DAModel
Definition: DAModel.H:59
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::DAObjFuncTotalPressure::areaSum_
scalar areaSum_
the area of all total pressure patches
Definition: DAObjFuncTotalPressure.H:39
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
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::DAObjFuncTotalPressure::DAObjFuncTotalPressure
DAObjFuncTotalPressure(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: DAObjFuncTotalPressure.C:19
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::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