DAResidualLaplacianFoam.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(DAResidualLaplacianFoam, 0);
16 addToRunTimeSelectionTable(DAResidual, DAResidualLaplacianFoam, 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  setResidualClassMemberScalar(T, dimensionSet(0, 0, -1, 1, 0, 0, 0))
28 
29 {
30  IOdictionary transportProperties(
31  IOobject(
32  "transportProperties",
33  mesh_.time().constant(),
34  mesh_,
35  IOobject::MUST_READ,
36  IOobject::NO_WRITE));
37 
38  dimensionedScalar DT("DT", dimViscosity, transportProperties);
39 
40  DT_ = DT.value();
41 }
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
45 {
46  /*
47  Description:
48  Clear all members to avoid memory leak because we will initalize
49  multiple objects of DAResidual. Here we need to delete all members
50  in the parent and child classes
51  */
52  TRes_.clear();
53 }
54 
55 void DAResidualLaplacianFoam::calcResiduals(const dictionary& options)
56 {
57  /*
58  Description:
59  This is the function to compute residuals.
60 
61  Input:
62  options.isPC: 1 means computing residuals for preconditioner matrix.
63  This essentially use the first order scheme for div(phi,U), div(phi,e)
64 
65  p_, T_, U_, phi_, etc: State variables in OpenFOAM
66 
67  Output:
68  URes_, pRes_, TRes_, phiRes_: residual field variables
69  */
70 
71  dimensionedScalar DT("DT", dimViscosity, DT_);
72 
73  fvScalarMatrix TEqn(
74  fvm::ddt(T_)
75  - fvm::laplacian(DT, T_));
76 
77  TRes_ = TEqn & T_;
78  normalizeResiduals(TRes);
79 }
80 
82 {
83  /*
84  Description:
85  Update the intermediate variables that depend on the state variables
86  */
87  // do nothing
88 }
89 
91 {
92  /*
93  Description:
94  Update the boundary condition for all the states in the selected solver
95  */
96 
97  T_.correctBoundaryConditions();
98 }
99 
100 } // End namespace Foam
101 
102 // ************************************************************************* //
Foam::DAResidualLaplacianFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualLaplacianFoam.C:90
Foam::DAResidualLaplacianFoam::clear
virtual void clear()
clear the members
Definition: DAResidualLaplacianFoam.C:44
TEqn
fvScalarMatrix TEqn(fvm::div(phi, T) - fvm::laplacian(alphaEff, T))
Foam::DAResidual::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAResidual.H:47
Foam::DAResidualLaplacianFoam::T_
volScalarField & T_
Definition: DAResidualLaplacianFoam.H:33
Foam::DAOption
Definition: DAOption.H:29
DT
dimensionedScalar & DT
Definition: createRefsLaplacian.H:6
DAResidualLaplacianFoam.H
daOption
DAOption daOption(mesh, pyOptions_)
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
Foam::DAResidualLaplacianFoam::TRes_
volScalarField TRes_
Definition: DAResidualLaplacianFoam.H:34
Foam::DAResidualLaplacianFoam::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAResidualLaplacianFoam.C:81
Foam::DAResidualLaplacianFoam::DAResidualLaplacianFoam
DAResidualLaplacianFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidualLaplacianFoam.C:19
Foam::DAResidualLaplacianFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualLaplacianFoam.C:55
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
Foam::DAResidual
Definition: DAResidual.H:35
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
Foam::DAResidualLaplacianFoam::DT_
scalar DT_
Definition: DAResidualLaplacianFoam.H:37
daModel
DAModel daModel(mesh, daOption)
daIndex
DAIndex daIndex(mesh, daOption, daModel)
setResidualClassMemberScalar
#define setResidualClassMemberScalar(stateName, stateUnit)
Definition: DAMacroFunctions.H:53