DAResidualTurboFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAResidualTurboFoam.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAResidualTurboFoam, 0);
16 addToRunTimeSelectionTable(DAResidual, DAResidualTurboFoam, dictionary);
17 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
18 
20  const word modelType,
21  const fvMesh& mesh,
22  const DAOption& daOption,
23  const DAModel& daModel,
24  const DAIndex& daIndex)
25  : DAResidual(modelType, mesh, daOption, daModel, daIndex),
26  // initialize and register state variables and their residuals, we use macros defined in macroFunctions.H
27  setResidualClassMemberVector(U, dimensionSet(1, -2, -2, 0, 0, 0, 0)),
28  setResidualClassMemberScalar(p, dimensionSet(1, -3, -1, 0, 0, 0, 0)),
29  setResidualClassMemberScalar(T, dimensionSet(1, -1, -3, 0, 0, 0, 0)),
31  thermo_(const_cast<fluidThermo&>(daModel.getThermo())),
32  URel_(const_cast<volVectorField&>(
33  mesh_.thisDb().lookupObject<volVectorField>("URel"))),
34  he_(thermo_.he()),
35  rho_(const_cast<volScalarField&>(
36  mesh_.thisDb().lookupObject<volScalarField>("rho"))),
37  alphat_(const_cast<volScalarField&>(
38  mesh_.thisDb().lookupObject<volScalarField>("alphat"))),
39  psi_(const_cast<volScalarField&>(
40  mesh_.thisDb().lookupObject<volScalarField>("thermo:psi"))),
41  daTurb_(const_cast<DATurbulenceModel&>(daModel.getDATurbulenceModel())),
42  // create simpleControl
43  simple_(const_cast<fvMesh&>(mesh)),
44  pressureControl_(p_, rho_, simple_.dict()),
45  MRF_(const_cast<IOMRFZoneListDF&>(
46  mesh_.thisDb().lookupObject<IOMRFZoneListDF>("MRFProperties")))
47 {
48 
49  // get molWeight and Cp from thermophysicalProperties
50  const IOdictionary& thermoDict = mesh.thisDb().lookupObject<IOdictionary>("thermophysicalProperties");
51  dictionary mixSubDict = thermoDict.subDict("mixture");
52  dictionary specieSubDict = mixSubDict.subDict("specie");
53  molWeight_ = specieSubDict.getScalar("molWeight");
54  dictionary thermodynamicsSubDict = mixSubDict.subDict("thermodynamics");
55  Cp_ = thermodynamicsSubDict.getScalar("Cp");
56 
57  if (daOption_.getOption<label>("debug"))
58  {
59  Info << "molWeight " << molWeight_ << endl;
60  Info << "Cp " << Cp_ << endl;
61  }
62 }
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
66 {
67  /*
68  Description:
69  Clear all members to avoid memory leak because we will initalize
70  multiple objects of DAResidual. Here we need to delete all members
71  in the parent and child classes
72  */
73  URes_.clear();
74  pRes_.clear();
75  TRes_.clear();
76  phiRes_.clear();
77 }
78 
79 void DAResidualTurboFoam::calcResiduals(const dictionary& options)
80 {
81  /*
82  Description:
83  This is the function to compute residuals.
84 
85  Input:
86  options.isPC: 1 means computing residuals for preconditioner matrix.
87  This essentially use the first order scheme for div(phi,U), div(phi,e)
88 
89  p_, T_, U_, phi_, etc: State variables in OpenFOAM
90 
91  Output:
92  URes_, pRes_, TRes_, phiRes_: residual field variables
93  */
94 
95  // We dont support MRF and fvOptions so all the related lines are commented
96  // out for now
97 
98  label isPC = options.getLabel("isPC");
99 
100  word divUScheme = "div(phi,U)";
101  word divHEScheme = "div(phi,e)";
102  word divPhidPScheme = "div(phid,p)";
103 
104  if (he_.name() == "h")
105  {
106  divHEScheme = "div(phi,h)";
107  }
108 
109  if (isPC)
110  {
111  divUScheme = "div(pc)";
112  divHEScheme = "div(pc)";
113  divPhidPScheme = "div(pc)";
114  }
115 
116  // ******** U Residuals **********
117  // copied and modified from UEqn.H
118 
119  tmp<fvVectorMatrix> tUEqn(
120  fvm::div(phi_, U_, divUScheme)
121  + MRF_.DDt(rho_, U_)
123  fvVectorMatrix& UEqn = tUEqn.ref();
124 
125  UEqn.relax();
126 
127  URes_ = (UEqn & U_) + fvc::grad(p_);
128  normalizeResiduals(URes);
129 
130  // ******** e Residuals **********
131  // copied and modified from EEqn.H
132  volSymmTensorField Teff = -daTurb_.devRhoReff();
133  volScalarField alphaEff("alphaEff", thermo_.alphaEff(alphat_));
134 
135  URel_ == U_;
137 
138  fvScalarMatrix EEqn(
139  fvm::div(phi_, he_, divHEScheme)
140  + (he_.name() == "e"
141  ? fvc::div(phi_, volScalarField("Ekp", 0.5 * magSqr(U_) + p_ / rho_))
142  : fvc::div(phi_, volScalarField("K", 0.5 * magSqr(U_))) - fvc::div(Teff.T() & U_) + fvc::div(p_ * (U_ - URel_)))
143  - fvm::Sp(fvc::div(phi_), he_)
144  - fvm::laplacian(alphaEff, he_));
145 
146  EEqn.relax();
147 
148  TRes_ = EEqn & he_;
149  normalizeResiduals(TRes);
150 
151  // ******** p and phi Residuals **********
152  // copied and modified from pEqn.H
153  volScalarField AU(UEqn.A());
154  volScalarField AtU(AU - UEqn.H1());
155  volVectorField HbyA("HbyA", U_);
156  HbyA = UEqn.H() / AU;
157 
158  volScalarField rAU(1.0 / UEqn.A());
159  tUEqn.clear();
160 
161  if (simple_.transonic())
162  {
163  surfaceScalarField phid(
164  "phid",
165  fvc::interpolate(psi_) * (fvc::interpolate(HbyA) & mesh_.Sf()));
166 
167  MRF_.makeRelative(fvc::interpolate(psi_), phid);
168 
169  fvScalarMatrix pEqn(
170  fvm::div(phid, p_, divPhidPScheme)
171  - fvm::laplacian(rho_ * rAU, p_));
172 
173  // for PC we do not include the div(phid, p) term, this improves the convergence
174  if (isPC && daOption_.getOption<label>("transonicPCOption") == 1)
175  {
176  pEqn -= fvm::div(phid, p_, divPhidPScheme);
177  }
178 
179  // Relax the pressure equation to maintain diagonal dominance
180  pEqn.relax();
181 
182  pEqn.setReference(pressureControl_.refCell(), pressureControl_.refValue());
183 
184  pRes_ = pEqn & p_;
185  normalizeResiduals(pRes);
186 
187  // ******** phi Residuals **********
188  // copied and modified from pEqn.H
189  if (isPC && daOption_.getOption<label>("transonicPCOption") == 2)
190  {
191  // transonic PC option 2, we ignore all the off-diagonal
192  // terms for the phiRes
193  phiRes_ = phi_;
194  }
195  else
196  {
197  phiRes_ == pEqn.flux() - phi_;
198  }
199 
200  // need to normalize phiRes
201  normalizePhiResiduals(phiRes);
202  }
203  else
204  {
205 
206  surfaceScalarField phiHbyA(
207  "phiHbyA",
208  fvc::interpolate(rho_ * HbyA) & mesh_.Sf());
209 
210  MRF_.makeRelative(fvc::interpolate(rho_), phiHbyA);
211 
212  adjustPhi(phiHbyA, U_, p_);
213  phiHbyA += fvc::interpolate(rho_ / AtU - rho_ / AU) * fvc::snGrad(p_) * mesh_.magSf();
214 
215  fvScalarMatrix pEqn(
216  fvc::div(phiHbyA)
217  - fvm::laplacian(rho_ / AtU, p_));
218 
219  pEqn.setReference(pressureControl_.refCell(), pressureControl_.refValue());
220 
221  pRes_ = pEqn & p_;
222  normalizeResiduals(pRes);
223 
224  // ******** phi Residuals **********
225  // copied and modified from pEqn.H
226  // TODO: the phiRes is not zero, need to fix
227  phiRes_ == phiHbyA + pEqn.flux() - phi_;
228  normalizePhiResiduals(phiRes);
229  }
230 }
231 
233 {
234  /*
235  Description:
236  Update the intermediate variables that depend on the state variables
237 
238  ********************** NOTE *****************
239  we assume hePsiThermo, pureMixture, perfectGas, hConst, and const transport
240  TODO: need to do this using built-in openfoam functions.
241 
242  we need to:
243  1, update psi based on T, psi=1/(R*T)
244  2, update rho based on p and psi, rho=psi*p
245  3, update E based on T, p and rho, E=Cp*T-p/rho
246  4, update velocity boundary based on MRF
247  */
248 
249  // 8314.4700665 gas constant in OpenFOAM
250  // src/OpenFOAM/global/constants/thermodynamic/thermodynamicConstants.H
251  scalar RR = Foam::constant::thermodynamic::RR;
252 
253  // R = RR/molWeight
254  // Foam::specie::R() function in src/thermophysicalModels/specie/specie/specieI.H
255  dimensionedScalar R(
256  "R1",
257  dimensionSet(0, 2, -2, -1, 0, 0, 0),
258  RR / molWeight_);
259 
260  // psi = 1/T/R
261  // see src/thermophysicalModels/specie/equationOfState/perfectGas/perfectGasI.H
262  psi_ = 1.0 / T_ / R;
263 
264  // rho = psi*p
265  // see src/thermophysicalModels/basic/psiThermo/psiThermo.C
266  rho_ = psi_ * p_;
267 
268  // **************** NOTE ****************
269  // need to relax rho to be consistent with the primal solver
270  // However, the rho.relax() will mess up perturbation
271  // That being said, we comment out the rho.relax() call to
272  // get the correct perturbed rho; however, the E residual will
273  // be a bit off compared with the ERes at the converged state
274  // from the primal solver. TODO. Need to figure out how to improve this
275  // **************** NOTE ****************
276  // rho_.relax();
277 
278  dimensionedScalar Cp(
279  "Cp1",
280  dimensionSet(0, 2, -2, -1, 0, 0, 0),
281  Cp_);
282 
283  // Hs = Cp*T
284  // see Hs() in src/thermophysicalModels/specie/thermo/hConst/hConstThermoI.H
285  // here the H departure EquationOfState::H(p, T) will be zero for perfectGas
286  // Es = Hs - p/rho = Hs - T * R;
287  // see Es() in src/thermophysicalModels/specie/thermo/thermo/thermoI.H
288  // **************** NOTE ****************
289  // See the comment from the rho.relax() call, if we write he_=Cp*T-p/rho, the
290  // accuracy of he_ may be impact by the inaccurate rho. So here we want to
291  // rewrite he_ as he_ = Cp * T_ - T_ * R instead, such that we dont include rho
292  // **************** NOTE ****************
293  if (he_.name() == "e")
294  {
295  he_ = Cp * T_ - T_ * R;
296  }
297  else
298  {
299  he_ = Cp * T_;
300  }
301  he_.correctBoundaryConditions();
302 
303  // NOTE: alphat is updated in the correctNut function in DATurbulenceModel child classes
304 
306 }
307 
309 {
310  /*
311  Description:
312  Update the boundary condition for all the states in the selected solver
313  */
314 
315  U_.correctBoundaryConditions();
316  p_.correctBoundaryConditions();
317  T_.correctBoundaryConditions();
318 
319 }
320 
321 } // End namespace Foam
322 
323 // ************************************************************************* //
phid
surfaceScalarField phid("phid",(fvc::interpolate(psi)/fvc::interpolate(rho)) *phiHbyA)
Foam::MRFZoneListDF::DDt
tmp< volVectorField > DDt(const volVectorField &U) const
Definition: MRFZoneListDF.C:120
Foam::DAResidualTurboFoam::pRes_
volScalarField pRes_
Definition: DAResidualTurboFoam.H:41
tUEqn
tmp< fvVectorMatrix > tUEqn(fvm::ddt(U)+fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
adjustPhi
adjustPhi(phiHbyA, U, p)
U
U
Definition: pEqnPimpleDyM.H:60
Foam::DAResidualTurboFoam::URes_
volVectorField URes_
Definition: DAResidualTurboFoam.H:38
Foam::DATurbulenceModel::divDevRhoReff
tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U)
divDev terms
Definition: DATurbulenceModel.C:324
Foam::DAResidual::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAResidual.H:47
Foam::DAResidualTurboFoam::clear
virtual void clear()
clear the members
Definition: DAResidualTurboFoam.C:65
Foam::DAResidualTurboFoam::TRes_
volScalarField TRes_
Definition: DAResidualTurboFoam.H:44
Foam::DAResidualTurboFoam::p_
volScalarField & p_
Definition: DAResidualTurboFoam.H:40
Foam::DAOption
Definition: DAOption.H:29
Foam::DAResidualTurboFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualTurboFoam.C:79
UEqn
fvVectorMatrix & UEqn
Definition: UEqnPimpleDyM.H:13
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAResidualTurboFoam::DAResidualTurboFoam
DAResidualTurboFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidualTurboFoam.C:19
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
Foam::DAResidualTurboFoam::simple_
simpleControl simple_
simpleControl object which will be initialized in this class
Definition: DAResidualTurboFoam.H:68
Foam::DAResidualTurboFoam::phi_
surfaceScalarField & phi_
Definition: DAResidualTurboFoam.H:46
Foam::DAOption::getOption
classType getOption(const word key) const
get an option from subDict and key
Definition: DAOption.H:92
HbyA
HbyA
Definition: pEqnRhoSimpleC.H:12
p
volScalarField & p
Definition: createRefsPimple.H:6
Foam::DAResidualTurboFoam::daTurb_
DATurbulenceModel & daTurb_
DATurbulenceModel object.
Definition: DAResidualTurboFoam.H:65
Foam::DAResidual::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAResidual.H:50
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
AU
volScalarField AU(UEqn.A())
setResidualClassMemberPhi
#define setResidualClassMemberPhi(stateName)
Definition: DAMacroFunctions.H:83
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
phiHbyA
phiHbyA
Definition: pEqnRhoSimpleC.H:27
Foam::DAIndex
Definition: DAIndex.H:32
DAResidualTurboFoam.H
setResidualClassMemberVector
#define setResidualClassMemberVector(stateName, stateUnit)
Definition: DAMacroFunctions.H:68
Foam::DAResidualTurboFoam::MRF_
IOMRFZoneListDF & MRF_
Multiple Reference Frame.
Definition: DAResidualTurboFoam.H:74
Foam::DAModel
Definition: DAModel.H:59
Foam::DAResidualTurboFoam::T_
volScalarField & T_
Definition: DAResidualTurboFoam.H:43
Foam
Definition: multiFreqScalarFvPatchField.C:144
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
Foam::DATurbulenceModel::devRhoReff
tmp< volSymmTensorField > devRhoReff() const
dev terms
Definition: DATurbulenceModel.C:306
AtU
volScalarField AtU(AU - UEqn.H1())
alphaEff
volScalarField alphaEff("alphaEff", turbulencePtr_->nu()/Pr+alphat)
Foam::IOMRFZoneListDF
Definition: IOMRFZoneListDF.H:49
Foam::DAResidualTurboFoam::alphat_
volScalarField & alphat_
Definition: DAResidualTurboFoam.H:60
Foam::MRFZoneListDF::makeRelative
void makeRelative(volVectorField &U) const
Definition: MRFZoneListDF.C:148
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
EEqn
fvScalarMatrix EEqn(fvm::div(phi, he)+(he.name()=="e" ? fvc::div(phi, volScalarField("Ekp", 0.5 *magSqr(U)+p/rho)) :fvc::div(phi, volScalarField("K", 0.5 *magSqr(U)))) - fvm::laplacian(turbulencePtr_->alphaEff(), he))
Foam::DAResidualTurboFoam::molWeight_
scalar molWeight_
Definition: DAResidualTurboFoam.H:78
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::DAResidualTurboFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualTurboFoam.C:308
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
he
volScalarField & he
Definition: EEqnRhoSimple.H:5
Foam::DAResidualTurboFoam::pressureControl_
pressureControl pressureControl_
pressureControl object to set ref pressure
Definition: DAResidualTurboFoam.H:71
rAU
volScalarField rAU(1.0/UEqn.A())
Foam::DAResidualTurboFoam::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAResidualTurboFoam.C:232
Foam::DAResidualTurboFoam::phiRes_
surfaceScalarField phiRes_
Definition: DAResidualTurboFoam.H:47
Foam::MRFZoneListDF::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Definition: MRFZoneListDF.C:207
daModel
DAModel daModel(mesh, daOption)
Foam::DAResidualTurboFoam::he_
volScalarField & he_
Definition: DAResidualTurboFoam.H:58
Foam::DAResidualTurboFoam::rho_
volScalarField & rho_
Definition: DAResidualTurboFoam.H:59
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAResidualTurboFoam::thermo_
fluidThermo & thermo_
thermophysical property
Definition: DAResidualTurboFoam.H:51
setResidualClassMemberScalar
#define setResidualClassMemberScalar(stateName, stateUnit)
Definition: DAMacroFunctions.H:53
Foam::DAResidualTurboFoam::U_
volVectorField & U_
Definition: DAResidualTurboFoam.H:37
Foam::DAResidualTurboFoam::Cp_
scalar Cp_
Definition: DAResidualTurboFoam.H:79
normalizePhiResiduals
#define normalizePhiResiduals(resName)
Definition: DAMacroFunctions.H:37
Foam::DAResidualTurboFoam::URel_
volVectorField & URel_
relative velocity
Definition: DAResidualTurboFoam.H:54
Teff
volSymmTensorField Teff
Definition: EEqnTurbo.H:3
Foam::DAResidualTurboFoam::psi_
volScalarField & psi_
Definition: DAResidualTurboFoam.H:61