DAPimpleFoam.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 DAPimpleFoam
8 
9  This class is modified from OpenFOAM's source code
10  applications/solvers/incompressible/pimpleFoam
11 
12  OpenFOAM: The Open Source CFD Toolbox
13 
14  Copyright (C): 2011-2016 OpenFOAM Foundation
15 
16  OpenFOAM License:
17 
18  OpenFOAM is free software: you can redistribute it and/or modify it
19  under the terms of the GNU General Public License as published by
20  the Free Software Foundation, either version 3 of the License, or
21  (at your option) any later version.
22 
23  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
24  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
25  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
26  for more details.
27 
28  You should have received a copy of the GNU General Public License
29  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #ifndef DAPimpleFoam_H
34 #define DAPimpleFoam_H
35 
36 #include "DASolver.H"
37 #include "addToRunTimeSelectionTable.H"
38 #include "singlePhaseTransportModel.H"
39 #include "turbulentTransportModel.H"
40 #include "pimpleControlDF.H"
41 #include "DAFvSource.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class DAPimpleFoam Declaration
50 \*---------------------------------------------------------------------------*/
51 
53  : public DASolver
54 {
55 
56 protected:
58  autoPtr<pimpleControlDF> pimplePtr_;
59 
61  autoPtr<volScalarField> pPtr_;
62 
64  autoPtr<volVectorField> UPtr_;
65 
67  autoPtr<surfaceScalarField> phiPtr_;
68 
70  autoPtr<singlePhaseTransportModel> laminarTransportPtr_;
71 
73  autoPtr<incompressible::turbulenceModel> turbulencePtr_;
74 
76  autoPtr<DATurbulenceModel> daTurbulenceModelPtr_;
77 
79  autoPtr<DAFvSource> daFvSourcePtr_;
80 
82  autoPtr<volVectorField> fvSourcePtr_;
83 
85  label hasFvSource_ = 0;
86 
88  scalar cumulativeContErr_ = 0.0;
89 
91  label pRefCell_ = 0;
92 
94  scalar pRefValue_ = 0.0;
95 
97  label IOInterval_ = 1;
98 
101 
103  autoPtr<dimensionedScalar> PrPtr_;
104 
106  autoPtr<dimensionedScalar> PrtPtr_;
107 
109  autoPtr<volScalarField> TPtr_;
110 
112  autoPtr<volScalarField> alphatPtr_;
113 
115  label hasTField_;
116 
117 public:
118  TypeName("DAPimpleFoam");
119  // Constructors
120 
121  //- Construct from components
122  DAPimpleFoam(
123  char* argsAll,
124  PyObject* pyOptions);
125 
126  //- Destructor
127  virtual ~DAPimpleFoam()
128  {
129  }
130 
132  virtual void initSolver();
133 
135  virtual label solvePrimal();
136 
138  virtual label solveAdjointFP(
139  Vec dFdW,
140  Vec psi);
141 
142  scalar calcAdjointResiduals(
143  const double* psi,
144  const double* dFdW,
145  double* adjRes);
146 };
147 
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149 
150 } // End namespace Foam
151 
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 
154 #endif
155 
156 // ************************************************************************* //
Foam::DAPimpleFoam::PrtPtr_
autoPtr< dimensionedScalar > PrtPtr_
Turbulence Prandtl pointer.
Definition: DAPimpleFoam.H:106
Foam::DAPimpleFoam::phiPtr_
autoPtr< surfaceScalarField > phiPtr_
surface flux field pointer
Definition: DAPimpleFoam.H:67
Foam::DAPimpleFoam::~DAPimpleFoam
virtual ~DAPimpleFoam()
Definition: DAPimpleFoam.H:127
Foam::DASolver
Definition: DASolver.H:54
Foam::DAPimpleFoam::pRefCell_
label pRefCell_
pressure referefence cell id
Definition: DAPimpleFoam.H:91
Foam::DAPimpleFoam::TPtr_
autoPtr< volScalarField > TPtr_
temperature field pointer
Definition: DAPimpleFoam.H:109
Foam::DAPimpleFoam::TypeName
TypeName("DAPimpleFoam")
DAFvSource.H
Foam::DAPimpleFoam::hasFvSource_
label hasFvSource_
whether to have fvSource term
Definition: DAPimpleFoam.H:85
Foam::DAPimpleFoam::fvSourcePtr_
autoPtr< volVectorField > fvSourcePtr_
fvSource term
Definition: DAPimpleFoam.H:82
Foam::DAPimpleFoam::DAPimpleFoam
DAPimpleFoam(char *argsAll, PyObject *pyOptions)
Definition: DAPimpleFoam.C:41
Foam::DAPimpleFoam::reduceIOWriteMesh_
label reduceIOWriteMesh_
whether to write mesh for the reduceIO
Definition: DAPimpleFoam.H:100
Foam::DAPimpleFoam::IOInterval_
label IOInterval_
the primal IO interval
Definition: DAPimpleFoam.H:97
Foam::DAPimpleFoam::solveAdjointFP
virtual label solveAdjointFP(Vec dFdW, Vec psi)
solve the adjoint equation using the fixed-point iteration method
Definition: DAPimpleFoam.C:268
Foam::DAPimpleFoam::cumulativeContErr_
scalar cumulativeContErr_
continuity error
Definition: DAPimpleFoam.H:88
pimpleControlDF.H
Foam::DAPimpleFoam::calcAdjointResiduals
scalar calcAdjointResiduals(const double *psi, const double *dFdW, double *adjRes)
Definition: DAPimpleFoam.C:242
Foam::DAPimpleFoam::PrPtr_
autoPtr< dimensionedScalar > PrPtr_
Prandtl number pointer.
Definition: DAPimpleFoam.H:103
Foam
Definition: checkGeometry.C:32
Foam::DAPimpleFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DAPimpleFoam.C:64
Foam::DAPimpleFoam::pimplePtr_
autoPtr< pimpleControlDF > pimplePtr_
pimple pointer
Definition: DAPimpleFoam.H:58
Foam::DAPimpleFoam::daTurbulenceModelPtr_
autoPtr< DATurbulenceModel > daTurbulenceModelPtr_
DATurbulenceModel pointer.
Definition: DAPimpleFoam.H:76
Foam::DAPimpleFoam::turbulencePtr_
autoPtr< incompressible::turbulenceModel > turbulencePtr_
turbulence pointer
Definition: DAPimpleFoam.H:73
Foam::DAPimpleFoam::pRefValue_
scalar pRefValue_
pressure reference value
Definition: DAPimpleFoam.H:94
Foam::DAPimpleFoam::daFvSourcePtr_
autoPtr< DAFvSource > daFvSourcePtr_
DASource pointer.
Definition: DAPimpleFoam.H:79
Foam::DAPimpleFoam::solvePrimal
virtual label solvePrimal()
solve the primal equations
Definition: DAPimpleFoam.C:119
DASolver.H
Foam::DAPimpleFoam::alphatPtr_
autoPtr< volScalarField > alphatPtr_
thermal diffusivity field pointer
Definition: DAPimpleFoam.H:112
Foam::DAPimpleFoam::laminarTransportPtr_
autoPtr< singlePhaseTransportModel > laminarTransportPtr_
laminar transport properties pointer
Definition: DAPimpleFoam.H:70
Foam::DAPimpleFoam
Definition: DAPimpleFoam.H:52
Foam::DAPimpleFoam::UPtr_
autoPtr< volVectorField > UPtr_
velocity field pointer
Definition: DAPimpleFoam.H:64
Foam::DAPimpleFoam::hasTField_
label hasTField_
whether to have the temperature field
Definition: DAPimpleFoam.H:115
Foam::DAPimpleFoam::pPtr_
autoPtr< volScalarField > pPtr_
pressure field pointer
Definition: DAPimpleFoam.H:61
psi
const volScalarField & psi
Definition: createRefsRhoPimple.H:14