DAResidualPimpleDyMFoam.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2 
3  DAFoam : Discrete Adjoint with OpenFOAM
4  Version : v4
5 
6  Description:
7  Child class for DAPimpleDyMFoam
8 
9 \*---------------------------------------------------------------------------*/
10 
11 #ifndef DAResidualPimpleDyMFoam_H
12 #define DAResidualPimpleDyMFoam_H
13 
14 #include "DAResidual.H"
15 #include "addToRunTimeSelectionTable.H"
16 #include "pimpleControl.H"
17 #include "adjustPhi.H"
18 
19 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
20 
21 namespace Foam
22 {
23 
24 /*---------------------------------------------------------------------------*\
25  Class DAResidualPimpleDyMFoam Declaration
26 \*---------------------------------------------------------------------------*/
27 
29  : public DAResidual
30 {
31 protected:
33 
34  volVectorField& U_;
35  volVectorField URes_;
36 
37  volScalarField& p_;
38  volScalarField pRes_;
39 
40  surfaceScalarField& phi_;
41  surfaceScalarField phiRes_;
42 
43  autoPtr<volScalarField> TResPtr_;
45 
47  volVectorField& fvSource_;
48 
51 
53  pimpleControl pimple_;
54 
56  label hasFvSource_ = 0;
57 
59  label hasTField_ = 0;
60 
61  scalar Pr_;
62 
63  scalar Prt_;
64 
65 public:
66  TypeName("DAPimpleDyMFoam");
67  // Constructors
68 
69  //- Construct from components
71  const word modelType,
72  const fvMesh& mesh,
73  const DAOption& daOption,
74  const DAModel& daModel,
75  const DAIndex& daIndex);
76 
77  //- Destructor
79  {
80  }
81 
82  // Members
83 
85  virtual void clear();
86 
88  virtual void calcResiduals(const dictionary& options);
89 
91  virtual void updateIntermediateVariables();
92 
94  virtual void correctBoundaryConditions();
95 
96  virtual void calcPCMatWithFvMatrix(Mat PCMat);
97 };
98 
99 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
100 
101 } // End namespace Foam
102 
103 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
104 
105 #endif
106 
107 // ************************************************************************* //
Foam::DAResidualPimpleDyMFoam
Definition: DAResidualPimpleDyMFoam.H:28
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
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::DAResidualPimpleDyMFoam::Pr_
scalar Pr_
Definition: DAResidualPimpleDyMFoam.H:61
Foam::DAResidualPimpleDyMFoam::~DAResidualPimpleDyMFoam
virtual ~DAResidualPimpleDyMFoam()
Definition: DAResidualPimpleDyMFoam.H:78
Foam::DAResidualPimpleDyMFoam::p_
volScalarField & p_
Definition: DAResidualPimpleDyMFoam.H:37
Foam::DAOption
Definition: DAOption.H:29
Foam::DAResidualPimpleDyMFoam::phi_
surfaceScalarField & phi_
Definition: DAResidualPimpleDyMFoam.H:40
Foam::DAResidualPimpleDyMFoam::U_
volVectorField & U_
Definition: DAResidualPimpleDyMFoam.H:34
Foam::DAResidualPimpleDyMFoam::pRes_
volScalarField pRes_
Definition: DAResidualPimpleDyMFoam.H:38
DAResidual.H
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
Foam::DAIndex
Definition: DAIndex.H:32
Foam::DAResidualPimpleDyMFoam::fvSource_
volVectorField & fvSource_
fvSource term
Definition: DAResidualPimpleDyMFoam.H:47
Foam::DAModel
Definition: DAModel.H:57
Foam
Definition: checkGeometry.C:32
Foam::DAResidualPimpleDyMFoam::URes_
volVectorField URes_
Definition: DAResidualPimpleDyMFoam.H:35
Foam::DAResidual
Definition: DAResidual.H:36
Foam::DATurbulenceModel
Definition: DATurbulenceModel.H:52
Foam::DAResidualPimpleDyMFoam::TypeName
TypeName("DAPimpleDyMFoam")
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
Foam::DAResidualPimpleDyMFoam::calcResiduals
virtual void calcResiduals(const dictionary &options)
compute residual
Definition: DAResidualPimpleDyMFoam.C:94
Foam::DAResidualPimpleDyMFoam::pimple_
pimpleControl pimple_
pimpleControl object which will be initialized in this class
Definition: DAResidualPimpleDyMFoam.H:53
Foam::DAResidualPimpleDyMFoam::correctBoundaryConditions
virtual void correctBoundaryConditions()
update the boundary condition for all the states in the selected solver
Definition: DAResidualPimpleDyMFoam.C:484