DAObjFuncMassFlowRate.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(DAObjFuncMassFlowRate, 0);
16 addToRunTimeSelectionTable(DAObjFunc, DAObjFuncMassFlowRate, 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 #ifdef CompressibleFlow
44  // setup the connectivity for total pressure, this is needed in Foam::DAJacCondFdW
45  word pName = "p";
46  if (mesh_.thisDb().foundObject<volScalarField>("p_rgh"))
47  {
48  pName = "p_rgh";
49  }
50  // for pressureInlet velocity, U depends on phi
51  objFuncConInfo_ = {{"U", "T", pName, "phi"}};
52 #endif
53 
54 #ifdef IncompressibleFlow
55  objFuncConInfo_ = {{"U"}};
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 mass flow rate M = integral( rho*U*dS )
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  // initialize faceValues to zero
89  forAll(objFuncFaceValues, idxI)
90  {
91  objFuncFaceValues[idxI] = 0.0;
92  }
93  // initialize objFunValue
94  objFuncValue = 0.0;
95 
96  const objectRegistry& db = mesh_.thisDb();
97  const volVectorField& U = db.lookupObject<volVectorField>("U");
98 
99  const volVectorField::Boundary& UBf = U.boundaryField();
100  const volScalarField::Boundary& rhoBf = rho_.boundaryField();
101 
102  forAll(objFuncFaceSources, idxI)
103  {
104  const label& objFuncFaceI = objFuncFaceSources[idxI];
105  label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
106  const label patchI = daIndex_.bFacePatchI[bFaceI];
107  const label faceI = daIndex_.bFaceFaceI[bFaceI];
108 
109  vector US = UBf[patchI][faceI];
110  vector Sf = mesh_.Sf().boundaryField()[patchI][faceI];
111  scalar rhoS = rhoBf[patchI][faceI];
112  scalar mfr = rhoS * (US & Sf);
113  objFuncFaceValues[idxI] = mfr * scale_;
114 
115  objFuncValue += objFuncFaceValues[idxI];
116  }
117 
118  // need to reduce the sum of force across all processors
119  reduce(objFuncValue, sumOp<scalar>());
120 
121  return;
122 }
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 } // End namespace Foam
127 
128 // ************************************************************************* //
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::DAObjFuncMassFlowRate::DAObjFuncMassFlowRate
DAObjFuncMassFlowRate(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: DAObjFuncMassFlowRate.C:19
DAObjFuncMassFlowRate.H
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
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAObjFuncMassFlowRate::rho_
volScalarField rho_
density
Definition: DAObjFuncMassFlowRate.H:36
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DAResidual
Definition: DAResidual.H:35
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
Foam::DAObjFuncMassFlowRate::calcObjFunc
virtual void calcObjFunc(const labelList &objFuncFaceSources, const labelList &objFuncCellSources, scalarList &objFuncFaceValues, scalarList &objFuncCellValues, scalar &objFuncValue)
calculate the value of objective function
Definition: DAObjFuncMassFlowRate.C:62
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::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