DAResidualRhoSimpleFoam.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(DAResidualRhoSimpleFoam, 0);
16 addToRunTimeSelectionTable(DAResidual, DAResidualRhoSimpleFoam, 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  fvSource_(const_cast<volVectorField&>(
32  mesh_.thisDb().lookupObject<volVectorField>("fvSource"))),
33  fvSourceEnergy_(const_cast<volScalarField&>(
34  mesh_.thisDb().lookupObject<volScalarField>("fvSourceEnergy"))),
35  fvOptions_(fv::options::New(mesh)),
36  thermo_(const_cast<fluidThermo&>(daModel.getThermo())),
37  he_(thermo_.he()),
38  rho_(const_cast<volScalarField&>(
39  mesh_.thisDb().lookupObject<volScalarField>("rho"))),
40  alphat_(const_cast<volScalarField&>(
41  mesh_.thisDb().lookupObject<volScalarField>("alphat"))),
42  psi_(const_cast<volScalarField&>(
43  mesh_.thisDb().lookupObject<volScalarField>("thermo:psi"))),
44  daTurb_(const_cast<DATurbulenceModel&>(daModel.getDATurbulenceModel())),
45  // create simpleControl
46  simple_(const_cast<fvMesh&>(mesh)),
47  pressureControl_(p_, rho_, simple_.dict()),
48  MRF_(const_cast<IOMRFZoneListDF&>(
49  mesh_.thisDb().lookupObject<IOMRFZoneListDF>("MRFProperties")))
50 {
51 
52  // initialize fvSource
53  const dictionary& allOptions = daOption.getAllOptions();
54  if (allOptions.subDict("fvSource").toc().size() != 0)
55  {
56  hasFvSource_ = 1;
57  }
58 
59  // get molWeight and Cp from thermophysicalProperties
60  const IOdictionary& thermoDict = mesh.thisDb().lookupObject<IOdictionary>("thermophysicalProperties");
61  dictionary mixSubDict = thermoDict.subDict("mixture");
62  dictionary specieSubDict = mixSubDict.subDict("specie");
63  molWeight_ = specieSubDict.getScalar("molWeight");
64  dictionary thermodynamicsSubDict = mixSubDict.subDict("thermodynamics");
65  Cp_ = thermodynamicsSubDict.getScalar("Cp");
66 
67  if (daOption_.getOption<label>("debug"))
68  {
69  Info << "molWeight " << molWeight_ << endl;
70  Info << "Cp " << Cp_ << endl;
71  }
72 
73  // this is just a dummy call because we need to run the constrain once
74  // to initialize fvOptions, before we can use it. Otherwise, we may get
75  // a seg fault when we call fvOptions_.correct(U_) in updateIntermediateVars
76  fvVectorMatrix UEqn(
77  fvm::div(phi_, U_)
78  - fvOptions_(rho_, U_));
79  fvOptions_.constrain(UEqn);
80 }
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
84 {
85  /*
86  Description:
87  Clear all members to avoid memory leak because we will initalize
88  multiple objects of DAResidual. Here we need to delete all members
89  in the parent and child classes
90  */
91  URes_.clear();
92  pRes_.clear();
93  TRes_.clear();
94  phiRes_.clear();
95 }
96 
97 void DAResidualRhoSimpleFoam::calcResiduals(const dictionary& options)
98 {
99  /*
100  Description:
101  This is the function to compute residuals.
102 
103  Input:
104  options.isPC: 1 means computing residuals for preconditioner matrix.
105  This essentially use the first order scheme for div(phi,U), div(phi,e)
106 
107  p_, T_, U_, phi_, etc: State variables in OpenFOAM
108 
109  Output:
110  URes_, pRes_, TRes_, phiRes_: residual field variables
111  */
112 
113  label isPC = options.getLabel("isPC");
114 
115  word divUScheme = "div(phi,U)";
116  word divHEScheme = "div(phi,e)";
117 
118  if (he_.name() == "h")
119  {
120  divHEScheme = "div(phi,h)";
121  }
122 
123  if (isPC)
124  {
125  divUScheme = "div(pc)";
126  divHEScheme = "div(pc)";
127  }
128 
129  // ******** U Residuals **********
130  // copied and modified from UEqn.H
131 
132  if (hasFvSource_)
133  {
134  DAFvSource& daFvSource(const_cast<DAFvSource&>(
135  mesh_.thisDb().lookupObject<DAFvSource>("DAFvSource")));
136  daFvSource.calcFvSource(fvSource_);
137  }
138 
139  tmp<fvVectorMatrix> tUEqn(
140  fvm::div(phi_, U_, divUScheme)
141  + MRF_.DDt(rho_, U_)
143  - fvSource_
144  - fvOptions_(rho_, U_));
145  fvVectorMatrix& UEqn = tUEqn.ref();
146 
147  UEqn.relax();
148 
149  fvOptions_.constrain(UEqn);
150 
151  URes_ = (UEqn & U_) + fvc::grad(p_);
152  normalizeResiduals(URes);
153 
154  // ******** e Residuals **********
155  // copied and modified from EEqn.H
156  volScalarField alphaEff("alphaEff", thermo_.alphaEff(alphat_));
157 
159 
160  fvScalarMatrix EEqn(
161  fvm::div(phi_, he_, divHEScheme)
162  + (he_.name() == "e"
163  ? fvc::div(phi_, volScalarField("Ekp", 0.5 * magSqr(U_) + p_ / rho_))
164  : fvc::div(phi_, volScalarField("K", 0.5 * magSqr(U_))))
165  - fvm::laplacian(alphaEff, he_)
167  - fvOptions_(rho_, he_));
168 
169  EEqn.relax();
170 
171  fvOptions_.constrain(EEqn);
172 
173  TRes_ = EEqn & he_;
174  normalizeResiduals(TRes);
175 
176  // ******** p and phi Residuals **********
177  // copied and modified from pEqn.H
178  volScalarField rAU(1.0 / UEqn.A());
179  surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho_ * rAU));
180  //volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
181  //***************** NOTE *******************
182  // constrainHbyA has been used since OpenFOAM-v1606; however, We do NOT use the constrainHbyA
183  // function in DAFoam because we found it significantly degrades the accuracy of shape derivatives.
184  // Basically, we should not constrain any variable because it will create discontinuity.
185  // Instead, we use the old implementation used in OpenFOAM-3.0+ and before
186  volVectorField HbyA("HbyA", U_);
187  HbyA = rAU * UEqn.H();
188  tUEqn.clear();
189 
190  surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(rho_) * fvc::flux(HbyA));
191 
192  MRF_.makeRelative(fvc::interpolate(rho_), phiHbyA);
193 
194  // Update the pressure BCs to ensure flux consistency
196 
197  // NOTE: we don't support transonic = true
198 
199  adjustPhi(phiHbyA, U_, p_);
200 
201  fvScalarMatrix pEqn(
202  fvc::div(phiHbyA)
203  - fvm::laplacian(rhorAUf, p_)
204  - fvOptions_(psi_, p_, rho_.name()));
205 
206  pEqn.setReference(pressureControl_.refCell(), pressureControl_.refValue());
207 
208  pRes_ = pEqn & p_;
209  normalizeResiduals(pRes);
210 
211  // ******** phi Residuals **********
212  // copied and modified from pEqn.H
213  phiRes_ = phiHbyA + pEqn.flux() - phi_;
214  normalizePhiResiduals(phiRes);
215 }
216 
218 {
219  /*
220  Description:
221  Update the intermediate variables that depend on the state variables
222 
223  ********************** NOTE *****************
224  we assume hePsiThermo, pureMixture, perfectGas, hConst, and const transport
225  TODO: need to do this using built-in openfoam functions.
226 
227  we need to:
228  1, update psi based on T, psi=1/(R*T)
229  2, update rho based on p and psi, rho=psi*p
230  3, update E based on T, p and rho, E=Cp*T-p/rho
231  4, update velocity boundary based on MRF
232  */
233 
234  // 8314.4700665 gas constant in OpenFOAM
235  // src/OpenFOAM/global/constants/thermodynamic/thermodynamicConstants.H
236  scalar RR = Foam::constant::thermodynamic::RR;
237 
238  // R = RR/molWeight
239  // Foam::specie::R() function in src/thermophysicalModels/specie/specie/specieI.H
240  dimensionedScalar R(
241  "R1",
242  dimensionSet(0, 2, -2, -1, 0, 0, 0),
243  RR / molWeight_);
244 
245  // psi = 1/T/R
246  // see src/thermophysicalModels/specie/equationOfState/perfectGas/perfectGasI.H
247  psi_ = 1.0 / T_ / R;
248 
249  // rho = psi*p
250  // see src/thermophysicalModels/basic/psiThermo/psiThermo.C
251  rho_ = psi_ * p_;
252 
253  // **************** NOTE ****************
254  // need to relax rho to be consistent with the primal solver
255  // However, the rho.relax() will mess up perturbation
256  // That being said, we comment out the rho.relax() call to
257  // get the correct perturbed rho; however, the E residual will
258  // be a bit off compared with the ERes at the converged state
259  // from the primal solver. TODO. Need to figure out how to improve this
260  // **************** NOTE ****************
261  // rho_.relax();
262 
263  dimensionedScalar Cp(
264  "Cp1",
265  dimensionSet(0, 2, -2, -1, 0, 0, 0),
266  Cp_);
267 
268  // Hs = Cp*T
269  // see Hs() in src/thermophysicalModels/specie/thermo/hConst/hConstThermoI.H
270  // here the H departure EquationOfState::H(p, T) will be zero for perfectGas
271  // Es = Hs - p/rho = Hs - T * R;
272  // see Es() in src/thermophysicalModels/specie/thermo/thermo/thermoI.H
273  // **************** NOTE ****************
274  // See the comment from the rho.relax() call, if we write he_=Cp*T-p/rho, the
275  // accuracy of he_ may be impact by the inaccurate rho. So here we want to
276  // rewrite he_ as he_ = Cp * T_ - T_ * R instead, such that we dont include rho
277  // **************** NOTE ****************
278  if (he_.name() == "e")
279  {
280  he_ = Cp * T_ - T_ * R;
281  }
282  else
283  {
284  he_ = Cp * T_;
285  }
286  he_.correctBoundaryConditions();
287 
288  // NOTE: alphat is updated in the correctNut function in DATurbulenceModel child classes
289 
291  fvOptions_.correct(U_);
292 }
293 
295 {
296  /*
297  Description:
298  Update the boundary condition for all the states in the selected solver
299  */
300 
301  U_.correctBoundaryConditions();
302  p_.correctBoundaryConditions();
303  T_.correctBoundaryConditions();
304 }
305 
306 } // End namespace Foam
307 
308 // ************************************************************************* //
Foam::DAFvSource
Definition: DAFvSource.H:34
Foam::MRFZoneListDF::DDt
tmp< volVectorField > DDt(const volVectorField &U) const
Definition: MRFZoneListDF.C:120
tUEqn
tmp< fvVectorMatrix > tUEqn(fvm::ddt(U)+fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
allOptions
const dictionary & allOptions
Definition: createRefsRhoSimpleC.H:15
Foam::DAResidualRhoSimpleFoam::hasFvSource_
label hasFvSource_
whether to have fvSource term
Definition: DAResidualRhoSimpleFoam.H:90
adjustPhi
adjustPhi(phiHbyA, U, p)
U
U
Definition: pEqnPimpleDyM.H:60
Foam::DAResidualRhoSimpleFoam::rho_
volScalarField & rho_
Definition: DAResidualRhoSimpleFoam.H:66
Foam::DAResidualRhoSimpleFoam::p_
volScalarField & p_
Definition: DAResidualRhoSimpleFoam.H:41
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::DAResidualRhoSimpleFoam::molWeight_
scalar molWeight_
Definition: DAResidualRhoSimpleFoam.H:85
Foam::DAResidualRhoSimpleFoam::Cp_
scalar Cp_
Definition: DAResidualRhoSimpleFoam.H:86
Foam::DAOption
Definition: DAOption.H:29
UEqn
fvVectorMatrix & UEqn
Definition: UEqnPimpleDyM.H:13
daOption
DAOption daOption(mesh, pyOptions_)
constrainPressure
constrainPressure(p, U, phiHbyA, rAtU(), MRF)
Foam::DAFvSource::calcFvSource
virtual void calcFvSource(volVectorField &fvSource)
compute the FvSource term
Definition: DAFvSource.C:90
Foam::DAResidualRhoSimpleFoam::phi_
surfaceScalarField & phi_
Definition: DAResidualRhoSimpleFoam.H:47
Foam::DAResidualRhoSimpleFoam::MRF_
IOMRFZoneListDF & MRF_
Multiple Reference Frame.
Definition: DAResidualRhoSimpleFoam.H:81
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
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::DAResidualRhoSimpleFoam::fvOptions_
fv::options & fvOptions_
fvOptions term
Definition: DAResidualRhoSimpleFoam.H:58
Foam::DAResidualRhoSimpleFoam::URes_
volVectorField URes_
Definition: DAResidualRhoSimpleFoam.H:39
Foam::DAResidual::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAResidual.H:50
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
setResidualClassMemberPhi
#define setResidualClassMemberPhi(stateName)
Definition: DAMacroFunctions.H:83
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
DAResidualRhoSimpleFoam.H
Foam::DAResidualRhoSimpleFoam::clear
virtual void clear()
clear the members
Definition: DAResidualRhoSimpleFoam.C:83
phiHbyA
phiHbyA
Definition: pEqnRhoSimpleC.H:27
Foam::DAResidualRhoSimpleFoam::psi_
volScalarField & psi_
Definition: DAResidualRhoSimpleFoam.H:68
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAResidualRhoSimpleFoam::thermo_
fluidThermo & thermo_
thermophysical property
Definition: DAResidualRhoSimpleFoam.H:61
setResidualClassMemberVector
#define setResidualClassMemberVector(stateName, stateUnit)
Definition: DAMacroFunctions.H:68
Foam::DAResidualRhoSimpleFoam::U_
volVectorField & U_
Definition: DAResidualRhoSimpleFoam.H:38
Foam::DAResidualRhoSimpleFoam::he_
volScalarField & he_
Definition: DAResidualRhoSimpleFoam.H:65
Foam::DAResidualRhoSimpleFoam::DAResidualRhoSimpleFoam
DAResidualRhoSimpleFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidualRhoSimpleFoam.C:19
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAResidualRhoSimpleFoam::pressureControl_
pressureControl pressureControl_
pressureControl object to set ref pressure
Definition: DAResidualRhoSimpleFoam.H:78
rhorAUf
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
alphaEff
volScalarField alphaEff("alphaEff", turbulencePtr_->nu()/Pr+alphat)
Foam::DAResidualRhoSimpleFoam::alphat_
volScalarField & alphat_
Definition: DAResidualRhoSimpleFoam.H:67
Foam::IOMRFZoneListDF
Definition: IOMRFZoneListDF.H:49
Foam::MRFZoneListDF::makeRelative
void makeRelative(volVectorField &U) const
Definition: MRFZoneListDF.C:148
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DAResidualRhoSimpleFoam::fvSourceEnergy_
volScalarField & fvSourceEnergy_
fvSource term for the energy equation
Definition: DAResidualRhoSimpleFoam.H:55
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::DAResidualRhoSimpleFoam::T_
volScalarField & T_
Definition: DAResidualRhoSimpleFoam.H:44
Foam::DAResidualRhoSimpleFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualRhoSimpleFoam.C:294
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
he
volScalarField & he
Definition: EEqnRhoSimple.H:5
rAU
volScalarField rAU(1.0/UEqn.A())
Foam::MRFZoneListDF::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Definition: MRFZoneListDF.C:207
daModel
DAModel daModel(mesh, daOption)
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAResidualRhoSimpleFoam::phiRes_
surfaceScalarField phiRes_
Definition: DAResidualRhoSimpleFoam.H:48
Foam::DAResidualRhoSimpleFoam::daTurb_
DATurbulenceModel & daTurb_
DATurbulenceModel object.
Definition: DAResidualRhoSimpleFoam.H:72
setResidualClassMemberScalar
#define setResidualClassMemberScalar(stateName, stateUnit)
Definition: DAMacroFunctions.H:53
Foam::DAResidualRhoSimpleFoam::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAResidualRhoSimpleFoam.C:217
Foam::DAResidualRhoSimpleFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualRhoSimpleFoam.C:97
normalizePhiResiduals
#define normalizePhiResiduals(resName)
Definition: DAMacroFunctions.H:37
Foam::DAResidualRhoSimpleFoam::fvSource_
volVectorField & fvSource_
fvSource term
Definition: DAResidualRhoSimpleFoam.H:52
Foam::DAResidualRhoSimpleFoam::TRes_
volScalarField TRes_
Definition: DAResidualRhoSimpleFoam.H:45
Foam::DAResidualRhoSimpleFoam::pRes_
volScalarField pRes_
Definition: DAResidualRhoSimpleFoam.H:42