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