DASimpleFoam.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 DASimpleFoam
8 
9  This class is modified from OpenFOAM's source code
10  applications/solvers/incompressible/simpleFoam
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 DASimpleFoam_H
34 #define DASimpleFoam_H
35 
36 #include "DASolver.H"
37 #include "addToRunTimeSelectionTable.H"
38 #include "singlePhaseTransportModel.H"
39 #include "turbulentTransportModel.H"
40 #include "simpleControl.H"
43 #include "DAFvSource.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class DASimpleFoam Declaration
52 \*---------------------------------------------------------------------------*/
53 
55  : public DASolver
56 {
57 
58 protected:
60  autoPtr<simpleControl> simplePtr_;
61 
63  autoPtr<volScalarField> pPtr_;
64 
66  autoPtr<volVectorField> UPtr_;
67 
69  autoPtr<surfaceScalarField> phiPtr_;
70 
72  autoPtr<volScalarField> alphaPorosityPtr_;
73 
75  autoPtr<singlePhaseTransportModel> laminarTransportPtr_;
76 
78  autoPtr<incompressible::turbulenceModel> turbulencePtr_;
79 
81  autoPtr<DATurbulenceModel> daTurbulenceModelPtr_;
82 
84  autoPtr<DAFvSource> daFvSourcePtr_;
85 
87  autoPtr<volVectorField> fvSourcePtr_;
88 
90  autoPtr<IOMRFZoneListDF> MRFPtr_;
91 
93  label hasFvSource_ = 0;
94 
96  scalar cumulativeContErr_ = 0.0;
97 
99  label pRefCell_ = 0;
100 
102  scalar pRefValue_ = 0.0;
103 
104 #ifdef CODI_AD_REVERSE
105  using Position = typename codi::RealReverse::Tape::Position;
107  Position adjResStart_;
108  Position adjResEnd_;
109  Position gradPStart_;
110  Position gradPEnd_;
111  Position pEqnFluxStart_;
112  Position pEqnFluxEnd_;
113  Position divPhiStart_;
114  Position divPhiEnd_;
115  Position fluxUStart_;
116  Position fluxUEnd_;
117 #endif
118 
120  scalar relaxUEqn_ = 1.0;
121  dictionary solverDictU_;
122  dictionary solverDictP_;
123 
124 public:
125  TypeName("DASimpleFoam");
126  // Constructors
127 
128  //- Construct from components
129  DASimpleFoam(
130  char* argsAll,
131  PyObject* pyOptions);
132 
133  //- Destructor
134  virtual ~DASimpleFoam()
135  {
136  }
137 
139  virtual void initSolver();
140 
142  virtual label solvePrimal(
143  const Vec xvVec,
144  Vec wVec);
145 
147  virtual label runFPAdj(
148  const Vec xvVec,
149  const Vec wVec,
150  Vec dFdW,
151  Vec psi);
152 
154  void invTranProdUEqn(
155  const volVectorField& mySource,
156  volVectorField& pseudoU);
157 
159  void invTranProdPEqn(
160  const volScalarField& mySource,
161  volScalarField& pseudoP);
162 
164  void calcLduResiduals(
165  volVectorField& URes,
166  volScalarField& pRes,
167  surfaceScalarField& phiRes);
168 
170  void vec2Fields(
171  const word mode,
172  Vec cVec,
173  volVectorField& UField,
174  volScalarField& pField,
175  surfaceScalarField& phiField,
176  volScalarField& nuTildaField);
177 
178  void calcAdjointResidual(
179  volVectorField& URes,
180  volScalarField& pRes,
181  surfaceScalarField& phiRes,
182  volScalarField& nuTildaRes,
183  Vec dFdW,
184  volVectorField& UPsi,
185  volScalarField& pPsi,
186  surfaceScalarField& phiPsi,
187  volScalarField& nuTildaPsi,
188  volVectorField& adjURes,
189  volScalarField& adjPRes,
190  surfaceScalarField& adjPhiRes,
191  volScalarField& adjNuTildaRes,
192  label& cnt);
193 
194  scalar L2norm(const volScalarField& v);
195 
196  vector L2norm(const volVectorField& U);
197 
198  scalar L2norm(const surfaceScalarField& Phi);
199 
200  void swap(List<scalar>& a, List<scalar>& b);
201 };
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 } // End namespace Foam
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 #endif
210 
211 // ************************************************************************* //
Foam::DASimpleFoam::phiPtr_
autoPtr< surfaceScalarField > phiPtr_
surface flux field pointer
Definition: DASimpleFoam.H:69
Foam::DASimpleFoam::DASimpleFoam
DASimpleFoam(char *argsAll, PyObject *pyOptions)
Definition: DASimpleFoam.C:41
Foam::DASimpleFoam::relaxUEqn_
scalar relaxUEqn_
fvSolution parameters for fixed-point adjoint
Definition: DASimpleFoam.H:120
Foam::DASimpleFoam::swap
void swap(List< scalar > &a, List< scalar > &b)
Definition: DASimpleFoam.C:1513
Foam::DASimpleFoam::pRefCell_
label pRefCell_
pressure referefence cell id
Definition: DASimpleFoam.H:99
U
U
Definition: pEqnPimpleDyM.H:60
DARegDbTurbulenceModelIncompressible.H
pseudoP
volScalarField pseudoP("pseudoP", p)
Foam::DASimpleFoam::invTranProdUEqn
void invTranProdUEqn(const volVectorField &mySource, volVectorField &pseudoU)
Inverse transpose product, MU^(-T)
Definition: DASimpleFoam.C:1049
Foam::DASimpleFoam::runFPAdj
virtual label runFPAdj(const Vec xvVec, const Vec wVec, Vec dFdW, Vec psi)
solve the adjoint equation using the fixed-point iteration method
Definition: DASimpleFoam.C:215
Foam::DASimpleFoam::fvSourcePtr_
autoPtr< volVectorField > fvSourcePtr_
fvSource term
Definition: DASimpleFoam.H:87
Foam::DASimpleFoam::solverDictU_
dictionary solverDictU_
Definition: DASimpleFoam.H:121
Foam::DASimpleFoam::daTurbulenceModelPtr_
autoPtr< DATurbulenceModel > daTurbulenceModelPtr_
DATurbulenceModel pointer.
Definition: DASimpleFoam.H:81
Foam::DASolver
Definition: DASolver.H:49
pseudoU
volVectorField pseudoU("pseudoU", U)
Foam::DASimpleFoam::alphaPorosityPtr_
autoPtr< volScalarField > alphaPorosityPtr_
alpha porosity
Definition: DASimpleFoam.H:72
DAFvSource.H
Foam::DASimpleFoam::laminarTransportPtr_
autoPtr< singlePhaseTransportModel > laminarTransportPtr_
laminar transport properties pointer
Definition: DASimpleFoam.H:75
Foam::DASimpleFoam::hasFvSource_
label hasFvSource_
whether to have fvSource term
Definition: DASimpleFoam.H:93
Foam::DASimpleFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DASimpleFoam.C:74
Foam::DASimpleFoam::TypeName
TypeName("DASimpleFoam")
Foam::DASimpleFoam::solverDictP_
dictionary solverDictP_
Definition: DASimpleFoam.H:122
Foam::DASimpleFoam::daFvSourcePtr_
autoPtr< DAFvSource > daFvSourcePtr_
DASource pointer.
Definition: DASimpleFoam.H:84
Foam::DASimpleFoam::cumulativeContErr_
scalar cumulativeContErr_
continuity error
Definition: DASimpleFoam.H:96
Foam::DASimpleFoam::~DASimpleFoam
virtual ~DASimpleFoam()
Definition: DASimpleFoam.H:134
Foam::DASimpleFoam::turbulencePtr_
autoPtr< incompressible::turbulenceModel > turbulencePtr_
turbulence pointer
Definition: DASimpleFoam.H:78
Foam::DASimpleFoam::L2norm
scalar L2norm(const volScalarField &v)
Definition: DASimpleFoam.C:1460
Foam
Definition: multiFreqScalarFvPatchField.C:144
Foam::DASimpleFoam::calcLduResiduals
void calcLduResiduals(volVectorField &URes, volScalarField &pRes, surfaceScalarField &phiRes)
calculate all the flow residuals except for the turbulence one
Definition: DASimpleFoam.C:1172
Foam::DASimpleFoam::pPtr_
autoPtr< volScalarField > pPtr_
pressure field pointer
Definition: DASimpleFoam.H:63
Foam::DASimpleFoam::invTranProdPEqn
void invTranProdPEqn(const volScalarField &mySource, volScalarField &pseudoP)
Inverse transpose product, Mp^(-T)
Definition: DASimpleFoam.C:1105
Foam::DASimpleFoam
Definition: DASimpleFoam.H:54
Foam::DASimpleFoam::simplePtr_
autoPtr< simpleControl > simplePtr_
simple pointer
Definition: DASimpleFoam.H:60
Foam::DASimpleFoam::vec2Fields
void vec2Fields(const word mode, Vec cVec, volVectorField &UField, volScalarField &pField, surfaceScalarField &phiField, volScalarField &nuTildaField)
assign coupled petsc vec to fields or vice versa
Definition: DASimpleFoam.C:942
Foam::DASimpleFoam::solvePrimal
virtual label solvePrimal(const Vec xvVec, Vec wVec)
solve the primal equations
Definition: DASimpleFoam.C:107
Foam::DASimpleFoam::calcAdjointResidual
void calcAdjointResidual(volVectorField &URes, volScalarField &pRes, surfaceScalarField &phiRes, volScalarField &nuTildaRes, Vec dFdW, volVectorField &UPsi, volScalarField &pPsi, surfaceScalarField &phiPsi, volScalarField &nuTildaPsi, volVectorField &adjURes, volScalarField &adjPRes, surfaceScalarField &adjPhiRes, volScalarField &adjNuTildaRes, label &cnt)
Definition: DASimpleFoam.C:1271
DASolver.H
Foam::DASimpleFoam::UPtr_
autoPtr< volVectorField > UPtr_
velocity field pointer
Definition: DASimpleFoam.H:66
Foam::DASimpleFoam::pRefValue_
scalar pRefValue_
pressure reference value
Definition: DASimpleFoam.H:102
Foam::DASimpleFoam::MRFPtr_
autoPtr< IOMRFZoneListDF > MRFPtr_
MRF pointer.
Definition: DASimpleFoam.H:90
psi
const volScalarField & psi
Definition: createRefsRhoSimpleC.H:14
dafoam_plot3dtransform.mode
mode
Definition: dafoam_plot3dtransform.py:21
DARegDbSinglePhaseTransportModel.H