DAResidualSimpleTFoam.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v3
5 
6  Description:
7  Child class for DASimpleTFoam
8 
9 \*---------------------------------------------------------------------------*/
10 
11 #ifndef DAResidualSimpleTFoam_H
12 #define DAResidualSimpleTFoam_H
13 
14 #include "DAResidual.H"
15 #include "addToRunTimeSelectionTable.H"
16 #include "simpleControl.H"
17 #include "adjustPhi.H"
18 #include "constrainPressure.H"
19 
20 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
21 
22 namespace Foam
23 {
24 
25 /*---------------------------------------------------------------------------*\
26  Class DAResidualSimpleTFoam Declaration
27 \*---------------------------------------------------------------------------*/
28 
30  : public DAResidual
31 {
32 protected:
34 
35  volVectorField& U_;
36  volVectorField URes_;
37 
38  volScalarField& p_;
39  volScalarField pRes_;
40 
41  volScalarField& T_;
42  volScalarField TRes_;
43 
44  surfaceScalarField& phi_;
45  surfaceScalarField phiRes_;
47 
49  volScalarField& alphaPorosity_;
50 
52  volVectorField& fvSource_;
53 
55  fv::options& fvOptions_;
56 
57  volScalarField& alphat_;
58 
59  scalar Pr_ = -9999.0;
60 
61  scalar Prt_ = -9999.0;
62 
65 
67  simpleControl simple_;
68 
71 
73  label hasFvSource_ = 0;
74 
75 public:
76  TypeName("DASimpleTFoam");
77  // Constructors
78 
79  //- Construct from components
81  const word modelType,
82  const fvMesh& mesh,
83  const DAOption& daOption,
84  const DAModel& daModel,
85  const DAIndex& daIndex);
86 
87  //- Destructor
89  {
90  }
91 
92  // Members
93 
95  virtual void clear();
96 
98  virtual void calcResiduals(const dictionary& options);
99 
101  virtual void updateIntermediateVariables();
102 
104  virtual void correctBoundaryConditions();
105 };
106 
107 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
108 
109 } // End namespace Foam
110 
111 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
112 
113 #endif
114 
115 // ************************************************************************* //
Foam::DAResidualSimpleTFoam::updateIntermediateVariables
virtual void updateIntermediateVariables()
update any intermediate variables that are dependent on state variables and are used in calcResiduals
Definition: DAResidualSimpleTFoam.C:205
Foam::DAResidualSimpleTFoam::alphaPorosity_
volScalarField & alphaPorosity_
alpha porosity term
Definition: DAResidualSimpleTFoam.H:49
Foam::DAResidualSimpleTFoam
Definition: DAResidualSimpleTFoam.H:29
Foam::DAResidualSimpleTFoam::clear
virtual void clear()
clear the members
Definition: DAResidualSimpleTFoam.C:74
Foam::DAResidualSimpleTFoam::simple_
simpleControl simple_
simpleControl object which will be initialized in this class
Definition: DAResidualSimpleTFoam.H:67
Foam::DAOption
Definition: DAOption.H:29
Foam::DAResidualSimpleTFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualSimpleTFoam.C:220
Foam::DAResidualSimpleTFoam::Pr_
scalar Pr_
Definition: DAResidualSimpleTFoam.H:59
Foam::DAResidualSimpleTFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualSimpleTFoam.C:88
daOption
DAOption daOption(mesh, pyOptions_)
Foam::DAResidualSimpleTFoam::daTurb_
DATurbulenceModel & daTurb_
DATurbulenceModel object.
Definition: DAResidualSimpleTFoam.H:64
Foam::DAResidualSimpleTFoam::U_
volVectorField & U_
Definition: DAResidualSimpleTFoam.H:35
Foam::DAResidualSimpleTFoam::DAResidualSimpleTFoam
DAResidualSimpleTFoam(const word modelType, const fvMesh &mesh, const DAOption &daOption, const DAModel &daModel, const DAIndex &daIndex)
Definition: DAResidualSimpleTFoam.C:19
Foam::DAResidualSimpleTFoam::phi_
surfaceScalarField & phi_
Definition: DAResidualSimpleTFoam.H:44
DAResidual.H
Foam::DAResidualSimpleTFoam::TypeName
TypeName("DASimpleTFoam")
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAResidualSimpleTFoam::MRF_
IOMRFZoneListDF & MRF_
Multiple Reference Frame.
Definition: DAResidualSimpleTFoam.H:70
Foam::DAModel
Definition: DAModel.H:59
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DAResidualSimpleTFoam::Prt_
scalar Prt_
Definition: DAResidualSimpleTFoam.H:61
Foam::DAResidualSimpleTFoam::hasFvSource_
label hasFvSource_
whether to have fvSource term
Definition: DAResidualSimpleTFoam.H:73
Foam::DAResidualSimpleTFoam::alphat_
volScalarField & alphat_
Definition: DAResidualSimpleTFoam.H:57
Foam::DAResidualSimpleTFoam::p_
volScalarField & p_
Definition: DAResidualSimpleTFoam.H:38
Foam::DAResidualSimpleTFoam::fvSource_
volVectorField & fvSource_
fvSource term
Definition: DAResidualSimpleTFoam.H:52
Foam::IOMRFZoneListDF
Definition: IOMRFZoneListDF.H:49
Foam::DAResidual
Definition: DAResidual.H:35
Foam::DAResidualSimpleTFoam::~DAResidualSimpleTFoam
virtual ~DAResidualSimpleTFoam()
Definition: DAResidualSimpleTFoam.H:88
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:61
Foam::DAResidualSimpleTFoam::phiRes_
surfaceScalarField phiRes_
Definition: DAResidualSimpleTFoam.H:45
daModel
DAModel daModel(mesh, daOption)
daIndex
DAIndex daIndex(mesh, daOption, daModel)
Foam::DAResidualSimpleTFoam::URes_
volVectorField URes_
Definition: DAResidualSimpleTFoam.H:36
Foam::DAResidualSimpleTFoam::pRes_
volScalarField pRes_
Definition: DAResidualSimpleTFoam.H:39
Foam::DAResidualSimpleTFoam::TRes_
volScalarField TRes_
Definition: DAResidualSimpleTFoam.H:42
Foam::DAResidualSimpleTFoam::fvOptions_
fv::options & fvOptions_
fvOptions term
Definition: DAResidualSimpleTFoam.H:55
Foam::DAResidualSimpleTFoam::T_
volScalarField & T_
Definition: DAResidualSimpleTFoam.H:41