DAResidualPimpleDyMFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6 \*---------------------------------------------------------------------------*/
7 
9 
10 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
11 
12 namespace Foam
13 {
14 
15 defineTypeNameAndDebug(DAResidualPimpleDyMFoam, 0);
16 addToRunTimeSelectionTable(DAResidual, DAResidualPimpleDyMFoam, 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  TResPtr_(nullptr),
31  fvSource_(const_cast<volVectorField&>(
32  mesh_.thisDb().lookupObject<volVectorField>("fvSource"))),
33  daTurb_(const_cast<DATurbulenceModel&>(daModel.getDATurbulenceModel())),
34  // create simpleControl
35  pimple_(const_cast<fvMesh&>(mesh))
36 {
37  // initialize fvSource
38  const dictionary& allOptions = daOption.getAllOptions();
39  if (allOptions.subDict("fvSource").toc().size() != 0)
40  {
41  hasFvSource_ = 1;
42  }
43 
44  // check whether to include the temperature field
45  hasTField_ = DAUtility::isFieldReadable(mesh, "T", "volScalarField");
46  if (hasTField_)
47  {
48 
49  TResPtr_.reset(new volScalarField(
50  IOobject(
51  "TRes",
52  mesh.time().timeName(),
53  mesh,
54  IOobject::NO_READ,
55  IOobject::NO_WRITE),
56  mesh,
57  dimensionedScalar("TRes", dimensionSet(0, 0, -1, 1, 0, 0, 0), 0.0),
58  zeroGradientFvPatchField<scalar>::typeName));
59 
60  // initialize the Prandtl number from transportProperties
61  IOdictionary transportProperties(
62  IOobject(
63  "transportProperties",
64  mesh.time().constant(),
65  mesh,
66  IOobject::MUST_READ,
67  IOobject::NO_WRITE,
68  false));
69  Pr_ = readScalar(transportProperties.lookup("Pr"));
70  Prt_ = readScalar(transportProperties.lookup("Prt"));
71  }
72 }
73 
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75 
77 {
78  /*
79  Description:
80  Clear all members to avoid memory leak because we will initalize
81  multiple objects of DAResidual. Here we need to delete all members
82  in the parent and child classes
83  */
84  URes_.clear();
85  pRes_.clear();
86  phiRes_.clear();
87 
88  if (hasTField_)
89  {
90  TResPtr_->clear();
91  }
92 }
93 
94 void DAResidualPimpleDyMFoam::calcResiduals(const dictionary& options)
95 {
96  /*
97  Description:
98  This is the function to compute residuals.
99 
100  Input:
101  options.isPC: 1 means computing residuals for preconditioner matrix.
102  This essentially use the first order scheme for div(phi,U)
103 
104  p_, U_, phi_, etc: State variables in OpenFOAM
105 
106  Output:
107  URes_, pRes_, phiRes_: residual field variables
108  */
109 
110  // We dont support MRF and fvOptions so all the related lines are commented
111  // out for now
112 
113  // ******** U Residuals **********
114  // copied and modified from UEqn.H
115 
116  word divUScheme = "div(phi,U)";
117 
118  label isPC = options.getLabel("isPC");
119 
120  if (isPC)
121  {
122  divUScheme = "div(pc)";
123  }
124 
125  if (hasFvSource_)
126  {
127  DAFvSource& daFvSource(const_cast<DAFvSource&>(
128  mesh_.thisDb().lookupObject<DAFvSource>("DAFvSource")));
129  // update the actuator source term
130  daFvSource.calcFvSource(fvSource_);
131  }
132 
133  fvVectorMatrix UEqn(
134  fvm::ddt(U_)
135  + fvm::div(phi_, U_, divUScheme)
137  - fvSource_);
138 
139  // NOTE: we need to call UEqn.relax here because it does some BC treatment, but we need to
140  // force the relaxation factor to be 1.0 because the last pimple loop does not use relaxation
141  UEqn.relax(1.0);
142 
143  URes_ = (UEqn & U_) + fvc::grad(p_);
144  normalizeResiduals(URes);
145 
146  // ******** p Residuals **********
147  // copied and modified from pEqn.H
148  // NOTE manually set pRefCell and pRefValue
149  label pRefCell = 0;
150  scalar pRefValue = 0.0;
151 
152  volScalarField rAU(1.0 / UEqn.A());
153  //***************** NOTE *******************
154  // constrainHbyA has been used since OpenFOAM-v1606; however, it may degrade the accuracy of derivatives
155  // because constraining variables will create discontinuity. Here we have a option to use the old
156  // implementation in OpenFOAM-3.0+ and before (no constraint for HbyA)
157  autoPtr<volVectorField> HbyAPtr = nullptr;
158  label useConstrainHbyA = daOption_.getOption<label>("useConstrainHbyA");
159  if (useConstrainHbyA)
160  {
161  HbyAPtr.reset(new volVectorField(constrainHbyA(rAU * UEqn.H(), U_, p_)));
162  }
163  else
164  {
165  HbyAPtr.reset(new volVectorField("HbyA", U_));
166  HbyAPtr() = rAU * UEqn.H();
167  }
168  volVectorField& HbyA = HbyAPtr();
169 
170  surfaceScalarField phiHbyA(
171  "phiHbyA",
172  fvc::flux(HbyA));
173 
174  tmp<volScalarField> rAtU(rAU);
175 
176  fvScalarMatrix pEqn(
177  fvm::laplacian(rAtU(), p_)
178  == fvc::div(phiHbyA));
179 
180  pEqn.setReference(pRefCell, pRefValue);
181 
182  pRes_ = pEqn & p_;
183  normalizeResiduals(pRes);
184 
185  // ******** phi Residuals **********
186  // copied and modified from pEqn.H
187  // NOTE: phiRes = phiHbyA - pEqn.flux() - phi; >> needs phi to be absolute instead of relative
188  // But we cannot simply << fvc::makeAbsolute(phi, U); >>, because UEqn still needs phi to be relative
189  // Therefore, we manually add the term fvc::meshPhi(U) to phiRes
190  phiRes_ = phiHbyA - pEqn.flux() - phi_ - fvc::meshPhi(U_);
191  // need to normalize phiRes
192  normalizePhiResiduals(phiRes);
193 
194  if (hasTField_)
195  {
196  volScalarField& alphat = const_cast<volScalarField&>(
197  mesh_.thisDb().lookupObject<volScalarField>("alphat"));
198 
199  volScalarField& T = const_cast<volScalarField&>(
200  mesh_.thisDb().lookupObject<volScalarField>("T"));
201 
202  volScalarField& TRes_ = TResPtr_();
203 
204  // ******** T Residuals **************
205  volScalarField alphaEff("alphaEff", daTurb_.nu() / Pr_ + alphat);
206 
207  fvScalarMatrix TEqn(
208  fvm::ddt(T)
209  + fvm::div(phi_, T)
210  - fvm::laplacian(alphaEff, T));
211 
212  TEqn.relax(1.0);
213 
214  TRes_ = TEqn & T;
215  normalizeResiduals(TRes);
216  }
217 }
218 
220 {
221  /*
222  Description:
223  Calculate the diagonal block of the preconditioner matrix dRdWTPC using the fvMatrix
224  */
225 
226  const labelUList& owner = mesh_.owner();
227  const labelUList& neighbour = mesh_.neighbour();
228 
229  PetscScalar val;
230 
231  dictionary normStateDict = daOption_.getAllOptions().subDict("normalizeStates");
232  wordList normResDict = daOption_.getOption<wordList>("normalizeResiduals");
233 
234  scalar UScaling = 1.0;
235  if (normStateDict.found("U"))
236  {
237  UScaling = normStateDict.getScalar("U");
238  }
239  scalar UResScaling = 1.0;
240 
241  fvVectorMatrix UEqn(
242  fvm::ddt(U_)
243  + fvm::div(phi_, U_, "div(pc)")
245  - fvSource_);
246 
247  UEqn.relax(1.0);
248 
249  // set diag
250  forAll(U_, cellI)
251  {
252  if (normResDict.found("URes"))
253  {
254  UResScaling = mesh_.V()[cellI];
255  }
256  for (label i = 0; i < 3; i++)
257  {
258  PetscInt rowI = daIndex_.getGlobalAdjointStateIndex("U", cellI, i);
259  PetscInt colI = rowI;
260  scalarField D = UEqn.D();
261  scalar val1 = D[cellI] * UScaling / UResScaling;
262  assignValueCheckAD(val, val1);
263  MatSetValues(PCMat, 1, &rowI, 1, &colI, &val, INSERT_VALUES);
264  }
265  }
266 
267  // set lower/owner
268  for (label faceI = 0; faceI < daIndex_.nLocalInternalFaces; faceI++)
269  {
270  label ownerCellI = owner[faceI];
271  label neighbourCellI = neighbour[faceI];
272 
273  if (normResDict.found("URes"))
274  {
275  UResScaling = mesh_.V()[neighbourCellI];
276  }
277 
278  for (label i = 0; i < 3; i++)
279  {
280  PetscInt rowI = daIndex_.getGlobalAdjointStateIndex("U", neighbourCellI, i);
281  PetscInt colI = daIndex_.getGlobalAdjointStateIndex("U", ownerCellI, i);
282  scalar val1 = UEqn.lower()[faceI] * UScaling / UResScaling;
283  assignValueCheckAD(val, val1);
284  MatSetValues(PCMat, 1, &colI, 1, &rowI, &val, INSERT_VALUES);
285  }
286  }
287 
288  // set upper/neighbour
289  for (label faceI = 0; faceI < daIndex_.nLocalInternalFaces; faceI++)
290  {
291  label ownerCellI = owner[faceI];
292  label neighbourCellI = neighbour[faceI];
293 
294  if (normResDict.found("URes"))
295  {
296  UResScaling = mesh_.V()[ownerCellI];
297  }
298 
299  for (label i = 0; i < 3; i++)
300  {
301  PetscInt rowI = daIndex_.getGlobalAdjointStateIndex("U", ownerCellI, i);
302  PetscInt colI = daIndex_.getGlobalAdjointStateIndex("U", neighbourCellI, i);
303  scalar val1 = UEqn.upper()[faceI] * UScaling / UResScaling;
304  assignValueCheckAD(val, val1);
305  MatSetValues(PCMat, 1, &colI, 1, &rowI, &val, INSERT_VALUES);
306  }
307  }
308 
309  label pRefCell = 0;
310  scalar pRefValue = 0.0;
311 
312  volScalarField rAU(1.0 / UEqn.A());
313  autoPtr<volVectorField> HbyAPtr = nullptr;
314  label useConstrainHbyA = daOption_.getOption<label>("useConstrainHbyA");
315  if (useConstrainHbyA)
316  {
317  HbyAPtr.reset(new volVectorField(constrainHbyA(rAU * UEqn.H(), U_, p_)));
318  }
319  else
320  {
321  HbyAPtr.reset(new volVectorField("HbyA", U_));
322  HbyAPtr() = rAU * UEqn.H();
323  }
324  volVectorField& HbyA = HbyAPtr();
325 
326  surfaceScalarField phiHbyA(
327  "phiHbyA",
328  fvc::flux(HbyA));
329 
330  tmp<volScalarField> rAtU(rAU);
331 
332  fvScalarMatrix pEqn(
333  fvm::laplacian(rAtU(), p_)
334  == fvc::div(phiHbyA));
335 
336  pEqn.setReference(pRefCell, pRefValue);
337 
338  // ********* p
339  scalar pScaling = 1.0;
340  if (normStateDict.found("p"))
341  {
342  pScaling = normStateDict.getScalar("p");
343  }
344  scalar pResScaling = 1.0;
345  // set diag
346  forAll(p_, cellI)
347  {
348  if (normResDict.found("pRes"))
349  {
350  pResScaling = mesh_.V()[cellI];
351  }
352 
353  PetscInt rowI = daIndex_.getGlobalAdjointStateIndex("p", cellI);
354  PetscInt colI = rowI;
355  scalarField D = pEqn.D();
356  scalar val1 = D[cellI] * pScaling / pResScaling;
357  assignValueCheckAD(val, val1);
358  MatSetValues(PCMat, 1, &rowI, 1, &colI, &val, INSERT_VALUES);
359  }
360 
361  // set lower/owner
362  for (label faceI = 0; faceI < daIndex_.nLocalInternalFaces; faceI++)
363  {
364  label ownerCellI = owner[faceI];
365  label neighbourCellI = neighbour[faceI];
366 
367  if (normResDict.found("pRes"))
368  {
369  pResScaling = mesh_.V()[neighbourCellI];
370  }
371 
372  PetscInt rowI = daIndex_.getGlobalAdjointStateIndex("p", neighbourCellI);
373  PetscInt colI = daIndex_.getGlobalAdjointStateIndex("p", ownerCellI);
374  scalar val1 = pEqn.lower()[faceI] * pScaling / pResScaling;
375  assignValueCheckAD(val, val1);
376  MatSetValues(PCMat, 1, &colI, 1, &rowI, &val, INSERT_VALUES);
377  }
378 
379  // set upper/neighbour
380  for (label faceI = 0; faceI < daIndex_.nLocalInternalFaces; faceI++)
381  {
382  label ownerCellI = owner[faceI];
383  label neighbourCellI = neighbour[faceI];
384 
385  if (normResDict.found("pRes"))
386  {
387  pResScaling = mesh_.V()[ownerCellI];
388  }
389 
390  PetscInt rowI = daIndex_.getGlobalAdjointStateIndex("p", ownerCellI);
391  PetscInt colI = daIndex_.getGlobalAdjointStateIndex("p", neighbourCellI);
392  scalar val1 = pEqn.upper()[faceI] * pScaling / pResScaling;
393  assignValueCheckAD(val, val1);
394  MatSetValues(PCMat, 1, &colI, 1, &rowI, &val, INSERT_VALUES);
395  }
396 
397  if (hasTField_)
398  {
399  volScalarField& alphat = const_cast<volScalarField&>(
400  mesh_.thisDb().lookupObject<volScalarField>("alphat"));
401 
402  volScalarField& T = const_cast<volScalarField&>(
403  mesh_.thisDb().lookupObject<volScalarField>("T"));
404 
405  volScalarField alphaEff("alphaEff", daTurb_.nu() / Pr_ + alphat);
406 
407  fvScalarMatrix TEqn(
408  fvm::ddt(T)
409  + fvm::div(phi_, T)
410  - fvm::laplacian(alphaEff, T));
411 
412  TEqn.relax(1.0);
413 
414  scalar TScaling = 1.0;
415  if (normStateDict.found("T"))
416  {
417  TScaling = normStateDict.getScalar("T");
418  }
419  scalar TResScaling = 1.0;
420  // set diag
421  forAll(T, cellI)
422  {
423  if (normResDict.found("TRes"))
424  {
425  TResScaling = mesh_.V()[cellI];
426  }
427 
428  PetscInt rowI = daIndex_.getGlobalAdjointStateIndex("T", cellI);
429  PetscInt colI = rowI;
430  scalarField D = TEqn.D();
431  scalar val1 = D[cellI] * TScaling / TResScaling;
432  assignValueCheckAD(val, val1);
433  MatSetValues(PCMat, 1, &rowI, 1, &colI, &val, INSERT_VALUES);
434  }
435 
436  // set lower/owner
437  for (label faceI = 0; faceI < daIndex_.nLocalInternalFaces; faceI++)
438  {
439  label ownerCellI = owner[faceI];
440  label neighbourCellI = neighbour[faceI];
441 
442  if (normResDict.found("TRes"))
443  {
444  TResScaling = mesh_.V()[neighbourCellI];
445  }
446 
447  PetscInt rowI = daIndex_.getGlobalAdjointStateIndex("T", neighbourCellI);
448  PetscInt colI = daIndex_.getGlobalAdjointStateIndex("T", ownerCellI);
449  scalar val1 = TEqn.lower()[faceI] * TScaling / TResScaling;
450  assignValueCheckAD(val, val1);
451  MatSetValues(PCMat, 1, &colI, 1, &rowI, &val, INSERT_VALUES);
452  }
453 
454  // set upper/neighbour
455  for (label faceI = 0; faceI < daIndex_.nLocalInternalFaces; faceI++)
456  {
457  label ownerCellI = owner[faceI];
458  label neighbourCellI = neighbour[faceI];
459 
460  if (normResDict.found("TRes"))
461  {
462  TResScaling = mesh_.V()[ownerCellI];
463  }
464 
465  PetscInt rowI = daIndex_.getGlobalAdjointStateIndex("T", ownerCellI);
466  PetscInt colI = daIndex_.getGlobalAdjointStateIndex("T", neighbourCellI);
467  scalar val1 = TEqn.upper()[faceI] * TScaling / TResScaling;
468  assignValueCheckAD(val, val1);
469  MatSetValues(PCMat, 1, &colI, 1, &rowI, &val, INSERT_VALUES);
470  }
471  }
472 }
473 
475 {
476  /*
477  Description:
478  Update the intermediate variables that depend on the state variables
479  */
480 
481  // nothing to update for DAPimpleDyMFoam
482 }
483 
485 {
486  /*
487  Description:
488  Update the boundary condition for all the states in the selected solver
489  */
490 
491  U_.correctBoundaryConditions();
492  p_.correctBoundaryConditions();
493  if (hasTField_)
494  {
495  volScalarField& T = const_cast<volScalarField&>(
496  mesh_.thisDb().lookupObject<volScalarField>("T"));
497  T.correctBoundaryConditions();
498  }
499 }
500 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
501 
502 } // End namespace Foam
503 
504 // ************************************************************************* //
Foam::DAFvSource
Definition: DAFvSource.H:34
U
volVectorField & U
Definition: createRefsPimpleDyM.H:7
pRefValue
scalar & pRefValue
Definition: createRefsPimpleDyM.H:12
Foam::DAResidualPimpleDyMFoam::Prt_
scalar Prt_
Definition: DAResidualPimpleDyMFoam.H:63
Foam::DAResidualPimpleDyMFoam::DAResidualPimpleDyMFoam
DAResidualPimpleDyMFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidualPimpleDyMFoam.C:19
Foam::DAResidualPimpleDyMFoam::calcPCMatWithFvMatrix
virtual void calcPCMatWithFvMatrix(Mat PCMat)
calculating the adjoint preconditioner matrix using fvMatrix
Definition: DAResidualPimpleDyMFoam.C:219
useConstrainHbyA
label useConstrainHbyA
Definition: pEqnPimpleDyM.H:7
Foam::DAResidualPimpleDyMFoam::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAResidualPimpleDyMFoam.C:474
Foam::DAResidual::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAResidual.H:49
Foam::DATurbulenceModel::divDevReff
tmp< fvVectorMatrix > divDevReff(volVectorField &U)
divDev terms
Definition: DATurbulenceModel.C:411
D
volVectorField & D
Definition: createRefsSolidDisplacement.H:8
Foam::DAResidualPimpleDyMFoam::Pr_
scalar Pr_
Definition: DAResidualPimpleDyMFoam.H:61
Foam::DAResidualPimpleDyMFoam::p_
volScalarField & p_
Definition: DAResidualPimpleDyMFoam.H:37
Foam::DAOption
Definition: DAOption.H:29
Foam::DAIndex::getGlobalAdjointStateIndex
label getGlobalAdjointStateIndex(const word stateName, const label idxI, const label comp=-1) const
get global adjoint index for a given state name, cell/face indxI and its component (optional,...
Definition: DAIndex.C:661
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
UEqn
fvVectorMatrix & UEqn
Definition: UEqnPimpleDyM.H:15
Foam::DAFvSource::calcFvSource
virtual void calcFvSource(volVectorField &fvSource)
compute the FvSource term
Definition: DAFvSource.C:90
Foam::DAResidualPimpleDyMFoam::phi_
surfaceScalarField & phi_
Definition: DAResidualPimpleDyMFoam.H:40
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
Foam::DAResidualPimpleDyMFoam::U_
volVectorField & U_
Definition: DAResidualPimpleDyMFoam.H:34
Foam::DAOption::getOption
classType getOption(const word key) const
get an option from subDict and key
Definition: DAOption.H:92
Foam::DAResidualPimpleDyMFoam::pRes_
volScalarField pRes_
Definition: DAResidualPimpleDyMFoam.H:38
Foam::DAResidual::daIndex_
const DAIndex & daIndex_
DAIndex.
Definition: DAResidual.H:58
pRefCell
label & pRefCell
Definition: createRefsPimpleDyM.H:11
Foam::DAOption::getAllOptions
const dictionary & getAllOptions() const
return a reference of allOptions_ dictionary
Definition: DAOption.H:56
Foam::DAResidual::daOption_
const DAOption & daOption_
DAOption object.
Definition: DAResidual.H:52
HbyA
volVectorField & HbyA
Definition: pEqnPimpleDyM.H:17
Foam::DAResidualPimpleDyMFoam::hasFvSource_
label hasFvSource_
whether to has fvSource term
Definition: DAResidualPimpleDyMFoam.H:56
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAResidualPimpleDyMFoam::phiRes_
surfaceScalarField phiRes_
Definition: DAResidualPimpleDyMFoam.H:41
setResidualClassMemberPhi
#define setResidualClassMemberPhi(stateName)
Definition: DAMacroFunctions.H:83
phiHbyA
phiHbyA
Definition: pEqnRhoSimpleC.H:36
Foam::DAIndex
Definition: DAIndex.H:32
p
volScalarField & p
Definition: createRefsPimpleDyM.H:6
assignValueCheckAD
#define assignValueCheckAD(valA, valB)
Definition: DAMacroFunctions.H:103
setResidualClassMemberVector
#define setResidualClassMemberVector(stateName, stateUnit)
Definition: DAMacroFunctions.H:68
Foam::DAResidualPimpleDyMFoam::fvSource_
volVectorField & fvSource_
fvSource term
Definition: DAResidualPimpleDyMFoam.H:47
Foam::DAModel
Definition: DAModel.H:57
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
phi
surfaceScalarField & phi
Definition: createRefsPimpleDyM.H:8
Foam
Definition: checkGeometry.C:32
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
rAtU
tmp< volScalarField > rAtU(rAU)
Foam::DAResidualPimpleDyMFoam::URes_
volVectorField URes_
Definition: DAResidualPimpleDyMFoam.H:35
allOptions
const dictionary & allOptions
Definition: createRefsPimpleDyM.H:14
Foam::DATurbulenceModel::nu
tmp< volScalarField > nu() const
get the nu field
Definition: DATurbulenceModel.C:267
alphaEff
volScalarField alphaEff("alphaEff", turbulencePtr_->nu()/Pr+alphat)
forAll
forAll(nuTilda1, cellI)
Definition: nuTilda1EqnIrkPimple.H:19
Foam::DAIndex::nLocalInternalFaces
label nLocalInternalFaces
local internal face size
Definition: DAIndex.H:98
Foam::DAResidual
Definition: DAResidual.H:36
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:52
alphat
volScalarField & alphat
Definition: TEqnPimpleDyM.H:5
rAU
volScalarField rAU(1.0/UEqn.A())
DAResidualPimpleDyMFoam.H
Foam::DAUtility::isFieldReadable
static label isFieldReadable(const fvMesh &mesh, const word fieldName, const word fieldType)
Definition: DAUtility.C:817
HbyAPtr
autoPtr< volVectorField > HbyAPtr
Definition: pEqnPimpleDyM.H:6
Foam::DAResidualPimpleDyMFoam::hasTField_
label hasTField_
whether to include the temperature field
Definition: DAResidualPimpleDyMFoam.H:59
Foam::DAResidualPimpleDyMFoam::daTurb_
DATurbulenceModel & daTurb_
DATurbulenceModel object.
Definition: DAResidualPimpleDyMFoam.H:50
Foam::DAResidualPimpleDyMFoam::clear
virtual void clear()
clear the members
Definition: DAResidualPimpleDyMFoam.C:76
Foam::DAResidualPimpleDyMFoam::TResPtr_
autoPtr< volScalarField > TResPtr_
Definition: DAResidualPimpleDyMFoam.H:43
setResidualClassMemberScalar
#define setResidualClassMemberScalar(stateName, stateUnit)
Definition: DAMacroFunctions.H:53
Foam::DAResidualPimpleDyMFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualPimpleDyMFoam.C:94
TEqn
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T))
normalizePhiResiduals
#define normalizePhiResiduals(resName)
Definition: DAMacroFunctions.H:37
Foam::DAResidualPimpleDyMFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualPimpleDyMFoam.C:484