DAResidualSimpleTFoam.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(DAResidualSimpleTFoam, 0);
16 addToRunTimeSelectionTable(DAResidual, DAResidualSimpleTFoam, 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(0, 1, -2, 0, 0, 0, 0)),
28  setResidualClassMemberScalar(p, dimensionSet(0, 0, -1, 0, 0, 0, 0)),
29  setResidualClassMemberScalar(T, dimensionSet(0, 0, -1, 1, 0, 0, 0)),
31  alphaPorosity_(const_cast<volScalarField&>(
32  mesh_.thisDb().lookupObject<volScalarField>("alphaPorosity"))),
33  fvSource_(const_cast<volVectorField&>(
34  mesh_.thisDb().lookupObject<volVectorField>("fvSource"))),
35  fvOptions_(fv::options::New(mesh)),
36  alphat_(const_cast<volScalarField&>(
37  mesh_.thisDb().lookupObject<volScalarField>("alphat"))),
38  daTurb_(const_cast<DATurbulenceModel&>(daModel.getDATurbulenceModel())),
39  // create simpleControl
40  simple_(const_cast<fvMesh&>(mesh)),
41  MRF_(const_cast<IOMRFZoneListDF&>(
42  mesh_.thisDb().lookupObject<IOMRFZoneListDF>("MRFProperties")))
43 {
44  // initialize fvSource
45  const dictionary& allOptions = daOption.getAllOptions();
46  if (allOptions.subDict("fvSource").toc().size() != 0)
47  {
48  hasFvSource_ = 1;
49  }
50 
51  // initialize the Prandtl number from transportProperties
52  IOdictionary transportProperties(
53  IOobject(
54  "transportProperties",
55  mesh.time().constant(),
56  mesh,
57  IOobject::MUST_READ,
58  IOobject::NO_WRITE,
59  false));
60  Pr_ = readScalar(transportProperties.lookup("Pr"));
61  Prt_ = readScalar(transportProperties.lookup("Prt"));
62 
63  // this is just a dummy call because we need to run the constrain once
64  // to initialize fvOptions, before we can use it. Otherwise, we may get
65  // a seg fault when we call fvOptions_.correct(U_) in updateIntermediateVars
66  fvVectorMatrix UEqn(
67  fvm::div(phi_, U_)
68  - fvOptions_(U_));
69  fvOptions_.constrain(UEqn);
70 }
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
75 {
76  /*
77  Description:
78  Clear all members to avoid memory leak because we will initalize
79  multiple objects of DAResidual. Here we need to delete all members
80  in the parent and child classes
81  */
82  URes_.clear();
83  pRes_.clear();
84  TRes_.clear();
85  phiRes_.clear();
86 }
87 
88 void DAResidualSimpleTFoam::calcResiduals(const dictionary& options)
89 {
90  /*
91  Description:
92  This is the function to compute residuals.
93 
94  Input:
95  options.isPC: 1 means computing residuals for preconditioner matrix.
96  This essentially use the first order scheme for div(phi,U)
97 
98  p_, U_, phi_, etc: State variables in OpenFOAM
99 
100  Output:
101  URes_, pRes_, phiRes_: residual field variables
102  */
103 
104  // We dont support MRF and fvOptions so all the related lines are commented
105  // out for now
106 
107  // ******** U Residuals **********
108  // copied and modified from UEqn.H
109 
110  word divUScheme = "div(phi,U)";
111 
112  label isPC = options.getLabel("isPC");
113 
114  if (isPC)
115  {
116  divUScheme = "div(pc)";
117  }
118 
119  if (hasFvSource_)
120  {
121  DAFvSource& daFvSource(const_cast<DAFvSource&>(
122  mesh_.thisDb().lookupObject<DAFvSource>("DAFvSource")));
123  daFvSource.calcFvSource(fvSource_);
124  }
125 
126  tmp<fvVectorMatrix> tUEqn(
127  fvm::div(phi_, U_, divUScheme)
128  + fvm::Sp(alphaPorosity_, U_)
129  + MRF_.DDt(U_)
131  - fvSource_
132  - fvOptions_(U_));
133  fvVectorMatrix& UEqn = tUEqn.ref();
134 
135  UEqn.relax();
136 
137  fvOptions_.constrain(UEqn);
138 
139  URes_ = (UEqn & U_) + fvc::grad(p_);
140  normalizeResiduals(URes);
141 
142  // ******** p Residuals **********
143  // copied and modified from pEqn.H
144  // NOTE manually set pRefCell and pRefValue
145  label pRefCell = 0;
146  scalar pRefValue = 0.0;
147 
148  volScalarField rAU(1.0 / UEqn.A());
149  //volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U_, p_));
150  //***************** NOTE *******************
151  // constrainHbyA has been used since OpenFOAM-v1606; however, We do NOT use the constrainHbyA
152  // function in DAFoam because we found it significantly degrades the accuracy of shape derivatives.
153  // Basically, we should not constrain any variable because it will create discontinuity.
154  // Instead, we use the old implementation used in OpenFOAM-3.0+ and before
155  volVectorField HbyA("HbyA", U_);
156  HbyA = rAU * UEqn.H();
157 
158  surfaceScalarField phiHbyA("phiHbyA", fvc::flux(HbyA));
159 
161 
162  adjustPhi(phiHbyA, U_, p_);
163 
164  tmp<volScalarField> rAtU(rAU);
165 
166  if (simple_.consistent())
167  {
168  rAtU = 1.0 / (1.0 / rAU - UEqn.H1());
169  phiHbyA += fvc::interpolate(rAtU() - rAU) * fvc::snGrad(p_) * mesh_.magSf();
170  HbyA -= (rAU - rAtU()) * fvc::grad(p_);
171  }
172 
173  tUEqn.clear();
174 
175  // Update the pressure BCs to ensure flux consistency
177 
178  fvScalarMatrix pEqn(
179  fvm::laplacian(rAtU(), p_)
180  == fvc::div(phiHbyA));
181  pEqn.setReference(pRefCell, pRefValue);
182 
183  pRes_ = pEqn & p_;
184  normalizeResiduals(pRes);
185 
186  // ******** phi Residuals **********
187  // copied and modified from pEqn.H
188  phiRes_ = phiHbyA - pEqn.flux() - phi_;
189  // need to normalize phiRes
190  normalizePhiResiduals(phiRes);
191 
192  // ******** T Residuals **************
193  volScalarField alphaEff("alphaEff", daTurb_.nu() / Pr_ + alphat_);
194 
195  fvScalarMatrix TEqn(
196  fvm::div(phi_, T_)
197  - fvm::laplacian(alphaEff, T_));
198 
199  TEqn.relax();
200 
201  TRes_ = TEqn & T_;
202  normalizeResiduals(TRes);
203 }
204 
206 {
207  /*
208  Description:
209  Update the intermediate variables that depend on the state variables
210  */
211 
212  alphat_ = daTurb_.getNut() / Prt_;
213  alphat_.correctBoundaryConditions();
214 
216 
217  fvOptions_.correct(U_);
218 }
219 
221 {
222  /*
223  Description:
224  Update the boundary condition for all the states in the selected solver
225  */
226 
227  U_.correctBoundaryConditions();
228  T_.correctBoundaryConditions();
229  p_.correctBoundaryConditions();
230 }
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace Foam
234 
235 // ************************************************************************* //
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
adjustPhi
adjustPhi(phiHbyA, U, p)
Foam::DAResidualSimpleTFoam::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAResidualSimpleTFoam.C:205
U
U
Definition: pEqnPimpleDyM.H:60
Foam::DAResidualSimpleTFoam::alphaPorosity_
volScalarField & alphaPorosity_
alpha porosity term
Definition: DAResidualSimpleTFoam.H:49
TEqn
fvScalarMatrix TEqn(fvm::div(phi, T) - fvm::laplacian(alphaEff, T))
Foam::DAResidual::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAResidual.H:47
Foam::DATurbulenceModel::divDevReff
tmp< fvVectorMatrix > divDevReff(volVectorField &U)
divDev terms
Definition: DATurbulenceModel.C:347
Foam::DAResidualSimpleTFoam::clear
virtual void clear()
clear the members
Definition: DAResidualSimpleTFoam.C:74
Foam::DAResidualSimpleTFoam::simple_
simpleControl simple_
simpleControl object which will be initialized in this class
Definition: DAResidualSimpleTFoam.H:67
pRefValue
scalar & pRefValue
Definition: createRefsPimple.H:12
Foam::DAOption
Definition: DAOption.H:29
Foam::DAResidualSimpleTFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualSimpleTFoam.C:220
Foam::DAResidualSimpleTFoam::Pr_
scalar Pr_
Definition: DAResidualSimpleTFoam.H:59
UEqn
fvVectorMatrix & UEqn
Definition: UEqnPimpleDyM.H:13
Foam::DAResidualSimpleTFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualSimpleTFoam.C:88
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::DAResidualSimpleTFoam::daTurb_
DATurbulenceModel & daTurb_
DATurbulenceModel object.
Definition: DAResidualSimpleTFoam.H:64
pRefCell
label & pRefCell
Definition: createRefsPimple.H:11
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
Foam::DAResidualSimpleTFoam::U_
volVectorField & U_
Definition: DAResidualSimpleTFoam.H:35
HbyA
HbyA
Definition: pEqnRhoSimpleC.H:12
Foam::DAResidualSimpleTFoam::DAResidualSimpleTFoam
DAResidualSimpleTFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidualSimpleTFoam.C:19
p
volScalarField & p
Definition: createRefsPimple.H:6
Foam::DAResidualSimpleTFoam::phi_
surfaceScalarField & phi_
Definition: DAResidualSimpleTFoam.H:44
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::DATurbulenceModel::getNut
tmp< volScalarField > getNut()
get the nut field
Definition: DATurbulenceModel.H:230
Foam::DAResidualSimpleTFoam::MRF_
IOMRFZoneListDF & MRF_
Multiple Reference Frame.
Definition: DAResidualSimpleTFoam.H:70
setResidualClassMemberVector
#define setResidualClassMemberVector(stateName, stateUnit)
Definition: DAMacroFunctions.H:68
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAResidualSimpleTFoam::Prt_
scalar Prt_
Definition: DAResidualSimpleTFoam.H:61
Foam::DAResidualSimpleTFoam::hasFvSource_
label hasFvSource_
whether to have fvSource term
Definition: DAResidualSimpleTFoam.H:73
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
Foam::DAResidualSimpleTFoam::alphat_
volScalarField & alphat_
Definition: DAResidualSimpleTFoam.H:57
Foam::DAResidualSimpleTFoam::p_
volScalarField & p_
Definition: DAResidualSimpleTFoam.H:38
alphaEff
volScalarField alphaEff("alphaEff", turbulencePtr_->nu()/Pr+alphat)
rAtU
tmp< volScalarField > rAtU(rAU)
Foam::DAResidualSimpleTFoam::fvSource_
volVectorField & fvSource_
fvSource term
Definition: DAResidualSimpleTFoam.H:52
Foam::IOMRFZoneListDF
Definition: IOMRFZoneListDF.H:49
Foam::MRFZoneListDF::makeRelative
void makeRelative(volVectorField &U) const
Definition: MRFZoneListDF.C:148
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:262
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
DAResidualSimpleTFoam.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DAResidualSimpleTFoam::phiRes_
surfaceScalarField phiRes_
Definition: DAResidualSimpleTFoam.H:45
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::DAResidualSimpleTFoam::URes_
volVectorField URes_
Definition: DAResidualSimpleTFoam.H:36
setResidualClassMemberScalar
#define setResidualClassMemberScalar(stateName, stateUnit)
Definition: DAMacroFunctions.H:53
Foam::DAResidualSimpleTFoam::pRes_
volScalarField pRes_
Definition: DAResidualSimpleTFoam.H:39
normalizePhiResiduals
#define normalizePhiResiduals(resName)
Definition: DAMacroFunctions.H:37
Foam::DAResidualSimpleTFoam::TRes_
volScalarField TRes_
Definition: DAResidualSimpleTFoam.H:42
Foam::DAResidualSimpleTFoam::fvOptions_
fv::options & fvOptions_
fvOptions term
Definition: DAResidualSimpleTFoam.H:55
Foam::DAResidualSimpleTFoam::T_
volScalarField & T_
Definition: DAResidualSimpleTFoam.H:41