DAFunctionTotalTemperatureRatio.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(DAFunctionTotalTemperatureRatio, 0);
16 addToRunTimeSelectionTable(DAFunction, DAFunctionTotalTemperatureRatio, 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 {
32 
33  functionDict_.readEntry<wordList>("inletPatches", inletPatches_);
34  functionDict_.readEntry<wordList>("outletPatches", outletPatches_);
35 
36  // get Cp from thermophysicalProperties
37  const IOdictionary& thermoDict = mesh_.thisDb().lookupObject<IOdictionary>("thermophysicalProperties");
38  dictionary mixSubDict = thermoDict.subDict("mixture");
39  dictionary thermodynamicsSubDict = mixSubDict.subDict("thermodynamics");
40  Cp_ = thermodynamicsSubDict.getScalar("Cp");
41  gamma_ = thermodynamicsSubDict.getScalar("gamma");
42 
43  if (daOption_.getOption<label>("debug"))
44  {
45  Info << "Cp " << Cp_ << endl;
46  Info << "gamma " << gamma_ << endl;
47  }
48 }
49 
52 {
53  /*
54  Description:
55  Calculate total temperature ratio, TTOut/TTIn
56  */
57 
58  // always calculate the area of all the inlet/outletPatches
59  areaSumInlet_ = 0.0;
60  areaSumOutlet_ = 0.0;
61  forAll(faceSources_, idxI)
62  {
63  const label& functionFaceI = faceSources_[idxI];
64  label bFaceI = functionFaceI - daIndex_.nLocalInternalFaces;
65  const label patchI = daIndex_.bFacePatchI[bFaceI];
66  const label faceI = daIndex_.bFaceFaceI[bFaceI];
67  word patchName = mesh_.boundaryMesh()[patchI].name();
68  if (inletPatches_.found(patchName))
69  {
70  areaSumInlet_ += mesh_.magSf().boundaryField()[patchI][faceI];
71  }
72  else if (outletPatches_.found(patchName))
73  {
74  areaSumOutlet_ += mesh_.magSf().boundaryField()[patchI][faceI];
75  }
76  else
77  {
78  FatalErrorIn(" ") << "inlet/outletPatches names are not in patches" << abort(FatalError);
79  }
80  }
81  reduce(areaSumInlet_, sumOp<scalar>());
82  reduce(areaSumOutlet_, sumOp<scalar>());
83 
84  // initialize objFunValue
85  scalar functionValue = 0.0;
86 
87  const objectRegistry& db = mesh_.thisDb();
88  const volScalarField& T = db.lookupObject<volScalarField>("T");
89  const volVectorField& U = db.lookupObject<volVectorField>("U");
90 
91  const scalar R = Cp_ - Cp_ / gamma_;
92 
93  const volScalarField::Boundary& TBf = T.boundaryField();
94  const volVectorField::Boundary& UBf = U.boundaryField();
95 
96  // we first compute the averaged inlet and outlet total pressure they will
97  // be used later for normalization
98  scalar TTIn = 0.0;
99  scalar TTOut = 0.0;
100  forAll(faceSources_, idxI)
101  {
102  const label& functionFaceI = faceSources_[idxI];
103  label bFaceI = functionFaceI - daIndex_.nLocalInternalFaces;
104  const label patchI = daIndex_.bFacePatchI[bFaceI];
105  const label faceI = daIndex_.bFaceFaceI[bFaceI];
106 
107  scalar TS = TBf[patchI][faceI];
108  scalar UMag = mag(UBf[patchI][faceI]);
109  scalar SfX = mesh_.magSf().boundaryField()[patchI][faceI];
110  scalar Ma2 = sqr(UMag / sqrt(gamma_ * R * TS));
111  scalar TT = TS * (1.0 + 0.5 * (gamma_ - 1.0) * Ma2);
112 
113  word patchName = mesh_.boundaryMesh()[patchI].name();
114  if (inletPatches_.found(patchName))
115  {
116  TTIn += TT * SfX / areaSumInlet_;
117  }
118  else if (outletPatches_.found(patchName))
119  {
120  TTOut += TT * SfX / areaSumOutlet_;
121  }
122  else
123  {
124  FatalErrorIn(" ") << "inlet/outletPatches names are not in patches" << abort(FatalError);
125  }
126  }
127  reduce(TTIn, sumOp<scalar>());
128  reduce(TTOut, sumOp<scalar>());
129 
130  functionValue = TTOut / TTIn;
131 
132  // check if we need to calculate refDiff.
133  this->calcRefVar(functionValue);
134 
135  return functionValue;
136 }
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 } // End namespace Foam
141 
142 // ************************************************************************* //
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::DAFunctionTotalTemperatureRatio::areaSumInlet_
scalar areaSumInlet_
the area of all inlet patches
Definition: DAFunctionTotalTemperatureRatio.H:32
Foam::DAOption
Definition: DAOption.H:29
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DAFunctionTotalTemperatureRatio::gamma_
scalar gamma_
Definition: DAFunctionTotalTemperatureRatio.H:45
Foam::DAOption::getOption
classType getOption(const word key) const
get an option from subDict and key
Definition: DAOption.H:92
DAFunctionTotalTemperatureRatio.H
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
Foam::DAFunctionTotalTemperatureRatio::calcFunction
virtual scalar calcFunction()
calculate the value of objective function
Definition: DAFunctionTotalTemperatureRatio.C:51
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAFunctionTotalTemperatureRatio::DAFunctionTotalTemperatureRatio
DAFunctionTotalTemperatureRatio(const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex, const word functionName)
Definition: DAFunctionTotalTemperatureRatio.C:19
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAModel
Definition: DAModel.H:57
Foam::DAFunction
Definition: DAFunction.H:31
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam::DAFunction::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAFunction.H:46
Foam
Definition: checkGeometry.C:32
Foam::DAFunctionTotalTemperatureRatio::inletPatches_
wordList inletPatches_
names of inlet patches
Definition: DAFunctionTotalTemperatureRatio.H:41
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DAFunctionTotalTemperatureRatio::Cp_
scalar Cp_
Definition: DAFunctionTotalTemperatureRatio.H:43
Foam::DAFunction::functionDict_
dictionary functionDict_
dictionary containing the information for the objective function
Definition: DAFunction.H:64
Foam::DAFunctionTotalTemperatureRatio::outletPatches_
wordList outletPatches_
names of outlet patches
Definition: DAFunctionTotalTemperatureRatio.H:38
Foam::DAFunctionTotalTemperatureRatio::areaSumOutlet_
scalar areaSumOutlet_
the area of all outlet patches
Definition: DAFunctionTotalTemperatureRatio.H:35
Foam::DAFunction::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAFunction.H:43
Foam::DAIndex::bFaceFaceI
labelList bFaceFaceI
Definition: DAIndex.H:116