DAResidualPimpleFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6 \*---------------------------------------------------------------------------*/
7 
8 #include "DAResidualPimpleFoam.H"
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAResidualPimpleFoam, 0);
16 addToRunTimeSelectionTable(DAResidual, DAResidualPimpleFoam, 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)),
30  fvSource_(const_cast<volVectorField&>(
31  mesh_.thisDb().lookupObject<volVectorField>("fvSource"))),
32  daTurb_(const_cast<DATurbulenceModel&>(daModel.getDATurbulenceModel())),
33  // create simpleControl
34  pimple_(const_cast<fvMesh&>(mesh))
35 {
36  // initialize fvSource
37  const dictionary& allOptions = daOption.getAllOptions();
38  if (allOptions.subDict("fvSource").toc().size() != 0)
39  {
40  hasFvSource_ = 1;
41  }
42 
43  mode_ = daOption.getSubDictOption<word>("unsteadyAdjoint", "mode");
44 }
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
49 {
50  /*
51  Description:
52  Clear all members to avoid memory leak because we will initalize
53  multiple objects of DAResidual. Here we need to delete all members
54  in the parent and child classes
55  */
56  URes_.clear();
57  pRes_.clear();
58  phiRes_.clear();
59 }
60 
61 void DAResidualPimpleFoam::calcResiduals(const dictionary& options)
62 {
63  /*
64  Description:
65  This is the function to compute residuals.
66 
67  Input:
68  options.isPC: 1 means computing residuals for preconditioner matrix.
69  This essentially use the first order scheme for div(phi,U)
70 
71  p_, U_, phi_, etc: State variables in OpenFOAM
72 
73  Output:
74  URes_, pRes_, phiRes_: residual field variables
75  */
76 
77  // We dont support MRF and fvOptions so all the related lines are commented
78  // out for now
79 
80  // ******** U Residuals **********
81  // copied and modified from UEqn.H
82 
83  word divUScheme = "div(phi,U)";
84 
85  label isPC = options.getLabel("isPC");
86 
87  if (isPC)
88  {
89  divUScheme = "div(pc)";
90  }
91 
92  if (hasFvSource_)
93  {
94  DAFvSource& daFvSource(const_cast<DAFvSource&>(
95  mesh_.thisDb().lookupObject<DAFvSource>("DAFvSource")));
96  // update the actuator source term
97  daFvSource.calcFvSource(fvSource_);
98  }
99 
100  fvVectorMatrix UEqn(
101  fvm::ddt(U_)
102  + fvm::div(phi_, U_, divUScheme)
104  - fvSource_);
105 
106  if (mode_ == "hybridAdjoint")
107  {
108  UEqn -= fvm::ddt(U_);
109  }
110 
111  UEqn.relax();
112 
113  URes_ = (UEqn & U_) + fvc::grad(p_);
114  normalizeResiduals(URes);
115 
116  // ******** p Residuals **********
117  // copied and modified from pEqn.H
118  // NOTE manually set pRefCell and pRefValue
119  label pRefCell = 0;
120  scalar pRefValue = 0.0;
121 
122  volScalarField rAU(1.0 / UEqn.A());
123  //volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U_, p_));
124  //***************** NOTE *******************
125  // constrainHbyA has been used since OpenFOAM-v1606; however, We do NOT use the constrainHbyA
126  // function in DAFoam because we found it significantly degrades the accuracy of shape derivatives.
127  // Basically, we should not constrain any variable because it will create discontinuity.
128  // Instead, we use the old implementation used in OpenFOAM-3.0+ and before
129  volVectorField HbyA("HbyA", U_);
130  HbyA = rAU * UEqn.H();
131 
132  surfaceScalarField phiHbyA(
133  "phiHbyA",
134  fvc::flux(HbyA)
135  + fvc::interpolate(rAU) * fvc::ddtCorr(U_, phi_));
136 
137  if (mode_ == "hybridAdjoint")
138  {
139  phiHbyA -= fvc::interpolate(rAU) * fvc::ddtCorr(U_, phi_);
140  }
141 
142  adjustPhi(phiHbyA, U_, p_);
143 
144  tmp<volScalarField> rAtU(rAU);
145 
146  if (pimple_.consistent())
147  {
148  rAtU = 1.0 / max(1.0 / rAU - UEqn.H1(), 0.1 / rAU);
149  phiHbyA +=
150  fvc::interpolate(rAtU() - rAU) * fvc::snGrad(p_) * mesh_.magSf();
151  HbyA -= (rAU - rAtU()) * fvc::grad(p_);
152  }
153 
154  fvScalarMatrix pEqn(
155  fvm::laplacian(rAtU(), p_)
156  == fvc::div(phiHbyA));
157 
158  pEqn.setReference(pRefCell, pRefValue);
159 
160  pRes_ = pEqn & p_;
161  normalizeResiduals(pRes);
162 
163  // ******** phi Residuals **********
164  // copied and modified from pEqn.H
165  phiRes_ = phiHbyA - pEqn.flux() - phi_;
166  // need to normalize phiRes
167  normalizePhiResiduals(phiRes);
168 }
169 
171 {
172  /*
173  Description:
174  Update the intermediate variables that depend on the state variables
175  */
176 
177  // nothing to update for DAPimpleFoam
178 }
179 
181 {
182  /*
183  Description:
184  Update the boundary condition for all the states in the selected solver
185  */
186 
187  U_.correctBoundaryConditions();
188  p_.correctBoundaryConditions();
189 }
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 } // End namespace Foam
193 
194 // ************************************************************************* //
DAResidualPimpleFoam.H
Foam::DAResidualPimpleFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualPimpleFoam.C:61
Foam::DAFvSource
Definition: DAFvSource.H:34
Foam::DAResidualPimpleFoam::daTurb_
DATurbulenceModel & daTurb_
DATurbulenceModel object.
Definition: DAResidualPimpleFoam.H:48
Foam::DAResidualPimpleFoam::p_
volScalarField & p_
Definition: DAResidualPimpleFoam.H:37
allOptions
const dictionary & allOptions
Definition: createRefsRhoSimpleC.H:15
adjustPhi
adjustPhi(phiHbyA, U, p)
U
U
Definition: pEqnPimpleDyM.H:60
Foam::DAResidualPimpleFoam::pimple_
pimpleControl pimple_
pimpleControl object which will be initialized in this class
Definition: DAResidualPimpleFoam.H:51
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
pRefValue
scalar & pRefValue
Definition: createRefsPimple.H:12
Foam::DAOption
Definition: DAOption.H:29
Foam::DAResidualPimpleFoam::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAResidualPimpleFoam.C:170
Foam::DAResidualPimpleFoam::mode_
word mode_
whether the hybrid adjoint or time accurate adjoint is active
Definition: DAResidualPimpleFoam.H:57
UEqn
fvVectorMatrix & UEqn
Definition: UEqnPimpleDyM.H:13
Foam::DAResidualPimpleFoam::fvSource_
volVectorField & fvSource_
fvSource term
Definition: DAResidualPimpleFoam.H:45
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAResidualPimpleFoam::U_
volVectorField & U_
Definition: DAResidualPimpleFoam.H:34
Foam::DAFvSource::calcFvSource
virtual void calcFvSource(volVectorField &fvSource)
compute the FvSource term
Definition: DAFvSource.C:90
Foam::DAResidualPimpleFoam::hasFvSource_
label hasFvSource_
whether to has fvSource term
Definition: DAResidualPimpleFoam.H:54
Foam::DAResidualPimpleFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualPimpleFoam.C:180
pRefCell
label & pRefCell
Definition: createRefsPimple.H:11
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
HbyA
HbyA
Definition: pEqnRhoSimpleC.H:12
p
volScalarField & p
Definition: createRefsPimple.H:6
Foam::DAResidualPimpleFoam::phi_
surfaceScalarField & phi_
Definition: DAResidualPimpleFoam.H:40
Foam::DAResidualPimpleFoam::URes_
volVectorField URes_
Definition: DAResidualPimpleFoam.H:35
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::DAResidualPimpleFoam::phiRes_
surfaceScalarField phiRes_
Definition: DAResidualPimpleFoam.H:41
Foam::DAResidualPimpleFoam::DAResidualPimpleFoam
DAResidualPimpleFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidualPimpleFoam.C:19
setResidualClassMemberVector
#define setResidualClassMemberVector(stateName, stateUnit)
Definition: DAMacroFunctions.H:68
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
rAtU
tmp< volScalarField > rAtU(rAU)
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
rAU
volScalarField rAU(1.0/UEqn.A())
daModel
DAModel daModel(mesh, daOption)
Foam::DAResidualPimpleFoam::clear
virtual void clear()
clear the members
Definition: DAResidualPimpleFoam.C:48
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAResidualPimpleFoam::pRes_
volScalarField pRes_
Definition: DAResidualPimpleFoam.H:38
setResidualClassMemberScalar
#define setResidualClassMemberScalar(stateName, stateUnit)
Definition: DAMacroFunctions.H:53
normalizePhiResiduals
#define normalizePhiResiduals(resName)
Definition: DAMacroFunctions.H:37