DAObjFuncTotalTemperatureRatio.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(DAObjFuncTotalTemperatureRatio, 0);
16 addToRunTimeSelectionTable(DAObjFunc, DAObjFuncTotalTemperatureRatio, 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  // setup the connectivity for total pressure, this is needed in Foam::DAJacCondFdW
42  // for pressureInlet velocity, U depends on phi
43  objFuncConInfo_ = {{"U", "T", "phi"}};
44 
45  objFuncDict_.readEntry<scalar>("scale", scale_);
46 
47  objFuncDict_.readEntry<wordList>("inletPatches", inletPatches_);
48  objFuncDict_.readEntry<wordList>("outletPatches", outletPatches_);
49 
50  // get Cp from thermophysicalProperties
51  const IOdictionary& thermoDict = mesh_.thisDb().lookupObject<IOdictionary>("thermophysicalProperties");
52  dictionary mixSubDict = thermoDict.subDict("mixture");
53  dictionary thermodynamicsSubDict = mixSubDict.subDict("thermodynamics");
54  Cp_ = thermodynamicsSubDict.getScalar("Cp");
55  gamma_ = thermodynamicsSubDict.getScalar("gamma");
56 
57  if (daOption_.getOption<label>("debug"))
58  {
59  Info << "Cp " << Cp_ << endl;
60  Info << "gamma " << gamma_ << endl;
61  }
62 }
63 
66  const labelList& objFuncFaceSources,
67  const labelList& objFuncCellSources,
68  scalarList& objFuncFaceValues,
69  scalarList& objFuncCellValues,
70  scalar& objFuncValue)
71 {
72  /*
73  Description:
74  Calculate total temperature ratio, TTIn/TTOut
75  NOTE: to enable coloring, we need to separate TTIn and TTOut for each face, we do:
76  d(TTInAvg/TTOutAvg)/dw = -TTOutAvg/TTInAvg^2 * dTTInAvg/dw + 1/TTInAvg * dTTOutAvg/dw
77  so here objDiscreteVal = -TTOutAvg/TTInAvg^2 * dTTIn_i*ds_i/dSIn for inlet patches
78  and objDiscreteVal = 1/TTInAvg * TTOut_i*ds_i/dSOut for outlet patches
79 
80  Input:
81  objFuncFaceSources: List of face source (index) for this objective
82 
83  objFuncCellSources: List of cell source (index) for this objective
84 
85  Output:
86  objFuncFaceValues: the discrete value of objective for each face source (index).
87  This will be used for computing df/dw in the adjoint.
88 
89  objFuncCellValues: the discrete value of objective on each cell source (index).
90  This will be used for computing df/dw in the adjoint.
91 
92  objFuncValue: the sum of objective, reduced across all processors and scaled by "scale"
93  */
94 
95  // calculate the area of all the inlet/outletPatches
96  if (areaSumInlet_ < 0.0 || areaSumOutlet_ < 0.0)
97  {
98  areaSumInlet_ = 0.0;
99  areaSumOutlet_ = 0.0;
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  word patchName = mesh_.boundaryMesh()[patchI].name();
107  if (DAUtility::isInList<word>(patchName, inletPatches_))
108  {
109  areaSumInlet_ += mesh_.magSf().boundaryField()[patchI][faceI];
110  }
111  else if (DAUtility::isInList<word>(patchName, outletPatches_))
112  {
113  areaSumOutlet_ += mesh_.magSf().boundaryField()[patchI][faceI];
114  }
115  else
116  {
117  FatalErrorIn(" ") << "inlet/outletPatches names are not in patches" << abort(FatalError);
118  }
119  }
120  reduce(areaSumInlet_, sumOp<scalar>());
121  reduce(areaSumOutlet_, sumOp<scalar>());
122  }
123 
124  // initialize faceValues to zero
125  forAll(objFuncFaceValues, idxI)
126  {
127  objFuncFaceValues[idxI] = 0.0;
128  }
129  // initialize objFunValue
130  objFuncValue = 0.0;
131 
132  const objectRegistry& db = mesh_.thisDb();
133  const volScalarField& T = db.lookupObject<volScalarField>("T");
134  const volVectorField& U = db.lookupObject<volVectorField>("U");
135 
136  const scalar R = Cp_ - Cp_ / gamma_;
137 
138  const volScalarField::Boundary& TBf = T.boundaryField();
139  const volVectorField::Boundary& UBf = U.boundaryField();
140 
141  // we first compute the averaged inlet and outlet total pressure they will
142  // be used later for normalization
143  scalar TTIn = 0.0;
144  scalar TTOut = 0.0;
145  forAll(objFuncFaceSources, idxI)
146  {
147  const label& objFuncFaceI = objFuncFaceSources[idxI];
148  label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
149  const label patchI = daIndex_.bFacePatchI[bFaceI];
150  const label faceI = daIndex_.bFaceFaceI[bFaceI];
151 
152  scalar TS = TBf[patchI][faceI];
153  scalar UMag = mag(UBf[patchI][faceI]);
154  scalar SfX = mesh_.magSf().boundaryField()[patchI][faceI];
155  scalar Ma2 = sqr(UMag / sqrt(gamma_ * R * TS));
156  scalar TT = TS * (1.0 + 0.5 * (gamma_ - 1.0) * Ma2);
157 
158  word patchName = mesh_.boundaryMesh()[patchI].name();
159  if (DAUtility::isInList<word>(patchName, inletPatches_))
160  {
161  TTIn += TT * SfX / areaSumInlet_;
162  }
163  else if (DAUtility::isInList<word>(patchName, outletPatches_))
164  {
165  TTOut += TT * SfX / areaSumOutlet_;
166  }
167  else
168  {
169  FatalErrorIn(" ") << "inlet/outletPatches names are not in patches" << abort(FatalError);
170  }
171  }
172  reduce(TTIn, sumOp<scalar>());
173  reduce(TTOut, sumOp<scalar>());
174 
175  // set Reference values
176  if (calcRefCoeffs)
177  {
178  TTInRef_ = TTIn;
179  TTOutRef_ = TTOut;
180  }
181 
182  objFuncValue = TTOut / TTIn;
183 
184  // pTotal ratio: Note we need special treatment to enable coloring
185  forAll(objFuncFaceSources, idxI)
186  {
187  const label& objFuncFaceI = objFuncFaceSources[idxI];
188  label bFaceI = objFuncFaceI - daIndex_.nLocalInternalFaces;
189  const label patchI = daIndex_.bFacePatchI[bFaceI];
190  const label faceI = daIndex_.bFaceFaceI[bFaceI];
191 
192  scalar TS = TBf[patchI][faceI];
193  scalar UMag = mag(UBf[patchI][faceI]);
194  scalar SfX = mesh_.magSf().boundaryField()[patchI][faceI];
195  scalar Ma2 = sqr(UMag / sqrt(gamma_ * R * TS));
196  scalar TT = TS * (1.0 + 0.5 * (gamma_ - 1.0) * Ma2);
197 
198  word patchName = mesh_.boundaryMesh()[patchI].name();
199  if (DAUtility::isInList<word>(patchName, inletPatches_))
200  {
201  objFuncFaceValues[idxI] = -TT * SfX / areaSumInlet_ * TTOutRef_ / TTInRef_ / TTInRef_;
202  }
203  else if (DAUtility::isInList<word>(patchName, outletPatches_))
204  {
205  objFuncFaceValues[idxI] = TT * SfX / areaSumOutlet_ / TTInRef_;
206  }
207  else
208  {
209  FatalErrorIn(" ") << "inlet/outletPatches names are not in patches" << abort(FatalError);
210  }
211  }
212 
213  return;
214 }
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace Foam
219 
220 // ************************************************************************* //
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
DAObjFuncTotalTemperatureRatio.H
Foam::DAIndex::bFacePatchI
labelList bFacePatchI
Definition: DAIndex.H:110
Foam::DAObjFuncTotalTemperatureRatio::outletPatches_
wordList outletPatches_
names of outlet patches
Definition: DAObjFuncTotalTemperatureRatio.H:39
Foam::DAOption
Definition: DAOption.H:29
Foam::DAObjFunc::daIndex_
const DAIndex & daIndex_
DAIndex object.
Definition: DAObjFunc.H:54
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAObjFuncTotalTemperatureRatio::calcObjFunc
virtual void calcObjFunc(const labelList &objFuncFaceSources, const labelList &objFuncCellSources, scalarList &objFuncFaceValues, scalarList &objFuncCellValues, scalar &objFuncValue)
calculate the value of objective function
Definition: DAObjFuncTotalTemperatureRatio.C:65
Foam::DAOption::getOption
classType getOption(const word key) const
get an option from subDict and key
Definition: DAOption.H:92
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::DAObjFunc::calcRefCoeffs
label calcRefCoeffs
whether to compute reference coefficients for special objFunc treatment such as totalPressureRatio
Definition: DAObjFunc.H:238
Foam::DAObjFuncTotalTemperatureRatio::TTOutRef_
scalar TTOutRef_
Definition: DAObjFuncTotalTemperatureRatio.H:46
Foam::DAObjFuncTotalTemperatureRatio::areaSumOutlet_
scalar areaSumOutlet_
the area of all outlet patches
Definition: DAObjFuncTotalTemperatureRatio.H:36
Foam::DAModel
Definition: DAModel.H:59
Foam::DAObjFuncTotalTemperatureRatio::TTInRef_
scalar TTInRef_
Definition: DAObjFuncTotalTemperatureRatio.H:44
Foam
Definition: multiFreqScalarFvPatchField.C:144
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
Foam::DAObjFuncTotalTemperatureRatio::inletPatches_
wordList inletPatches_
names of inlet patches
Definition: DAObjFuncTotalTemperatureRatio.H:42
Foam::DAObjFunc::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAObjFunc.H:48
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DAObjFuncTotalTemperatureRatio::Cp_
scalar Cp_
Definition: DAObjFuncTotalTemperatureRatio.H:48
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::DAObjFuncTotalTemperatureRatio::DAObjFuncTotalTemperatureRatio
DAObjFuncTotalTemperatureRatio(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: DAObjFuncTotalTemperatureRatio.C:19
Foam::DAObjFunc::objFuncConInfo_
List< List< word > > objFuncConInfo_
the connectivity information for the objective function
Definition: DAObjFunc.H:93
Foam::DAObjFuncTotalTemperatureRatio::areaSumInlet_
scalar areaSumInlet_
the area of all inlet patches
Definition: DAObjFuncTotalTemperatureRatio.H:33
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAObjFuncTotalTemperatureRatio::gamma_
scalar gamma_
Definition: DAObjFuncTotalTemperatureRatio.H:50
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