DAResidualHeatTransferFoam.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(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 solidProperties(
33  IOobject(
34  "solidProperties",
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  solidProperties));
45 
46  const dictionary& allOptions = daOption.getAllOptions();
47  if (allOptions.subDict("fvSource").toc().size() != 0)
48  {
49  hasFvSource_ = 1;
50  }
51 }
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
55 {
56  /*
57  Description:
58  Clear all members to avoid memory leak because we will initalize
59  multiple objects of DAResidual. Here we need to delete all members
60  in the parent and child classes
61  */
62  TRes_.clear();
63 }
64 
65 void DAResidualHeatTransferFoam::calcResiduals(const dictionary& options)
66 {
67  /*
68  Description:
69  This is the function to compute residuals.
70 
71  Input:
72  options.isPC: 1 means computing residuals for preconditioner matrix.
73  This essentially use the first order scheme for div(phi,U), div(phi,e)
74 
75  p_, T_, U_, phi_, etc: State variables in OpenFOAM
76 
77  Output:
78  URes_, pRes_, TRes_, phiRes_: residual field variables
79  */
80 
81  if (hasFvSource_)
82  {
83  DAFvSource& daFvSource(const_cast<DAFvSource&>(
84  mesh_.thisDb().lookupObject<DAFvSource>("DAFvSource")));
85  daFvSource.calcFvSource(fvSource_);
86  }
87 
88  fvScalarMatrix TEqn(
89  fvm::laplacian(kPtr_(), T_)
90  + fvSource_);
91 
92  TRes_ = TEqn & T_;
93  normalizeResiduals(TRes);
94 }
95 
97 {
98  /*
99  Description:
100  Update the intermediate variables that depend on the state variables
101  */
102  // do nothing
103 }
104 
106 {
107  /*
108  Description:
109  Update the boundary condition for all the states in the selected solver
110  */
111 
112  T_.correctBoundaryConditions();
113 }
114 
115 } // End namespace Foam
116 
117 // ************************************************************************* //
Foam::DAResidualHeatTransferFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualHeatTransferFoam.C:105
Foam::DAFvSource
Definition: DAFvSource.H:34
Foam::DAResidualHeatTransferFoam::DAResidualHeatTransferFoam
DAResidualHeatTransferFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidualHeatTransferFoam.C:19
Foam::DAResidual::mesh_
const fvMesh & mesh_
fvMesh
Definition: DAResidual.H:49
DAResidualHeatTransferFoam.H
Foam::DAResidualHeatTransferFoam::hasFvSource_
label hasFvSource_
Definition: DAResidualHeatTransferFoam.H:40
Foam::DAOption
Definition: DAOption.H:29
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(DAFunction, DAFunctionForce, dictionary)
Foam::DAFvSource::calcFvSource
virtual void calcFvSource(volVectorField &fvSource)
compute the FvSource term
Definition: DAFvSource.C:90
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:96
Foam::DAOption::getAllOptions
const dictionary & getAllOptions() const
return a reference of allOptions_ dictionary
Definition: DAOption.H:56
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAResidualHeatTransferFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualHeatTransferFoam.C:65
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAModel
Definition: DAModel.H:57
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(DAFunction, 0)
Foam
Definition: checkGeometry.C:32
T
volScalarField & T
Definition: createRefsHeatTransfer.H:5
allOptions
const dictionary & allOptions
Definition: createRefsPimpleDyM.H:14
Foam::DAResidual
Definition: DAResidual.H:36
Foam::DAResidualHeatTransferFoam::TRes_
volScalarField TRes_
Definition: DAResidualHeatTransferFoam.H:34
Foam::DAResidualHeatTransferFoam::clear
virtual void clear()
clear the members
Definition: DAResidualHeatTransferFoam.C:54
Foam::DAResidualHeatTransferFoam::fvSource_
volScalarField & fvSource_
Definition: DAResidualHeatTransferFoam.H:35
setResidualClassMemberScalar
#define setResidualClassMemberScalar(stateName, stateUnit)
Definition: DAMacroFunctions.H:53
TEqn
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T))
Foam::DAResidualHeatTransferFoam::T_
volScalarField & T_
Definition: DAResidualHeatTransferFoam.H:33
Foam::DAResidualHeatTransferFoam::kPtr_
autoPtr< dimensionedScalar > kPtr_
Definition: DAResidualHeatTransferFoam.H:38