DAResidualHeatTransferFoam.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(DAResidualHeatTransferFoam, 0);
16 addToRunTimeSelectionTable(DAResidual, DAResidualHeatTransferFoam, 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, dimPower / dimLength / dimLength / dimLength),
28  fvSource_(const_cast<volScalarField&>(mesh_.thisDb().lookupObject<volScalarField>("fvSource"))),
29  kPtr_(nullptr)
30 
31 {
32  IOdictionary transportProperties(
33  IOobject(
34  "transportProperties",
35  mesh_.time().constant(),
36  mesh_,
37  IOobject::MUST_READ,
38  IOobject::NO_WRITE));
39 
40  kPtr_.reset(
41  new dimensionedScalar(
42  "k",
43  dimPower / dimLength / dimTemperature,
44  transportProperties));
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  TRes_.clear();
57 }
58 
59 void DAResidualHeatTransferFoam::calcResiduals(const dictionary& options)
60 {
61  /*
62  Description:
63  This is the function to compute residuals.
64 
65  Input:
66  options.isPC: 1 means computing residuals for preconditioner matrix.
67  This essentially use the first order scheme for div(phi,U), div(phi,e)
68 
69  p_, T_, U_, phi_, etc: State variables in OpenFOAM
70 
71  Output:
72  URes_, pRes_, TRes_, phiRes_: residual field variables
73  */
74 
75  fvScalarMatrix TEqn(
76  fvm::laplacian(kPtr_(), T_)
77  + fvSource_);
78 
79  TRes_ = TEqn & T_;
80  normalizeResiduals(TRes);
81 }
82 
84 {
85  /*
86  Description:
87  Update the intermediate variables that depend on the state variables
88  */
89  // do nothing
90 }
91 
93 {
94  /*
95  Description:
96  Update the boundary condition for all the states in the selected solver
97  */
98 
99  T_.correctBoundaryConditions();
100 }
101 
102 } // End namespace Foam
103 
104 // ************************************************************************* //
Foam::DAResidualHeatTransferFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualHeatTransferFoam.C:92
Foam::DAResidualHeatTransferFoam::DAResidualHeatTransferFoam
DAResidualHeatTransferFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidualHeatTransferFoam.C:19
TEqn
fvScalarMatrix TEqn(fvm::div(phi, T) - fvm::laplacian(alphaEff, T))
Foam::DAResidual::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAResidual.H:47
DAResidualHeatTransferFoam.H
Foam::DAOption
Definition: DAOption.H:29
daOption
DAOption daOption(mesh, pyOptions_)
normalizeResiduals
#define normalizeResiduals(resName)
Definition: DAMacroFunctions.H:28
Foam::DAResidualHeatTransferFoam::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAResidualHeatTransferFoam.C:83
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAResidualHeatTransferFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualHeatTransferFoam.C:59
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::DAResidualHeatTransferFoam::TRes_
volScalarField TRes_
Definition: DAResidualHeatTransferFoam.H:34
Foam::DAResidualHeatTransferFoam::clear
virtual void clear()
clear the members
Definition: DAResidualHeatTransferFoam.C:48
Foam::DAResidualHeatTransferFoam::fvSource_
volScalarField & fvSource_
Definition: DAResidualHeatTransferFoam.H:35
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFvSource, 0)
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFvSource, DAFvSourceActuatorDisk, dictionary)
daModel
DAModel daModel(mesh, daOption)
daIndex
DAIndex daIndex(mesh, daOption, daModel)
setResidualClassMemberScalar
#define setResidualClassMemberScalar(stateName, stateUnit)
Definition: DAMacroFunctions.H:53
Foam::DAResidualHeatTransferFoam::T_
volScalarField & T_
Definition: DAResidualHeatTransferFoam.H:33
Foam::DAResidualHeatTransferFoam::kPtr_
autoPtr< dimensionedScalar > kPtr_
Definition: DAResidualHeatTransferFoam.H:38