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