DAFunctionMassFlowRate.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6 \*---------------------------------------------------------------------------*/
7 
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAFunctionMassFlowRate, 0);
16 addToRunTimeSelectionTable(DAFunction, DAFunctionMassFlowRate, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  const fvMesh& mesh,
21  const DAOption& daOption,
22  const DAModel& daModel,
23  const DAIndex& daIndex,
24  const word functionName)
25  : DAFunction(
26  mesh,
27  daOption,
28  daModel,
29  daIndex,
30  functionName),
31  daTurb_(const_cast<DATurbulenceModel&>(daModel.getDATurbulenceModel()))
32 {
33 }
34 
37 {
38  /*
39  Description:
40  Calculate mass flow rate M = integral( rho*U*dS )
41  */
42 
43  // initialize objFunValue
44  scalar functionValue = 0.0;
45 
46  const objectRegistry& db = mesh_.thisDb();
47  const volVectorField& U = db.lookupObject<volVectorField>("U");
48 
49  const volVectorField::Boundary& UBf = U.boundaryField();
50  volScalarField rho = daTurb_.rho();
51  const volScalarField::Boundary& rhoBf = rho.boundaryField();
52 
53  forAll(faceSources_, idxI)
54  {
55  const label& functionFaceI = faceSources_[idxI];
56  label bFaceI = functionFaceI - daIndex_.nLocalInternalFaces;
57  const label patchI = daIndex_.bFacePatchI[bFaceI];
58  const label faceI = daIndex_.bFaceFaceI[bFaceI];
59 
60  vector US = UBf[patchI][faceI];
61  vector Sf = mesh_.Sf().boundaryField()[patchI][faceI];
62  scalar rhoS = rhoBf[patchI][faceI];
63  scalar mfr = rhoS * (US & Sf) * scale_;
64 
65  functionValue += mfr;
66  }
67 
68  // need to reduce the sum of force across all processors
69  reduce(functionValue, sumOp<scalar>());
70 
71  // check if we need to calculate refDiff.
72  this->calcRefVar(functionValue);
73 
74  return functionValue;
75 }
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 } // End namespace Foam
80 
81 // ************************************************************************* //
DAFunctionMassFlowRate.H
U
volVectorField & U
Definition: createRefsPimpleDyM.H:7
Foam::DAFunction::faceSources_
labelList faceSources_
a sorted list of all face sources for the objective function
Definition: DAFunction.H:67
Foam::DAIndex::bFacePatchI
labelList bFacePatchI
Definition: DAIndex.H:110
Foam::DAOption
Definition: DAOption.H:29
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DAFunctionMassFlowRate::calcFunction
virtual scalar calcFunction()
calculate the value of objective function
Definition: DAFunctionMassFlowRate.C:36
Foam::DAFunctionMassFlowRate::DAFunctionMassFlowRate
DAFunctionMassFlowRate(const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const word functionName)
Definition: DAFunctionMassFlowRate.C:19
Foam::DATurbulenceModel::rho
tmp< volScalarField > rho() const
get the density field
Definition: DATurbulenceModel.C:294
Foam::DAFunction::calcRefVar
void calcRefVar(scalar &functionValue)
calculate (var-ref)^2
Definition: DAFunction.C:216
Foam::DAFunction::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAFunction.H:52
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAModel
Definition: DAModel.H:57
Foam::DAFunction
Definition: DAFunction.H:31
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam
Definition: checkGeometry.C:32
Foam::DAFunctionMassFlowRate::daTurb_
DATurbulenceModel & daTurb_
DATurbulenceModel object.
Definition: DAFunctionMassFlowRate.H:32
rho
volScalarField & rho
Definition: createRefsRhoPimple.H:8
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:52
Foam::DAFunction::scale_
scalar scale_
scale of the objective function
Definition: DAFunction.H:73
Foam::DAFunction::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAFunction.H:43
Foam::DAIndex::bFaceFaceI
labelList bFaceFaceI
Definition: DAIndex.H:116