DAIrkPimpleFoam.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 DAIrkPimpleFoam
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 DAIrkPimpleFoam_H
34 #define DAIrkPimpleFoam_H
35 
36 #include "DASolver.H"
37 #include "addToRunTimeSelectionTable.H"
38 #include "singlePhaseTransportModel.H"
39 #include "turbulentTransportModel.H"
40 //#include "pimpleControlDF.H" // A modified version of pimpleControl, Basically, we disable the output to the screen
41 #include "pimpleControl.H" // origional pimpleControl
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class DAIrkPimpleFoam Declaration
50 \*---------------------------------------------------------------------------*/
51 
53  : public DASolver
54 {
55 
56 protected:
57  // Radau23 coefficients
58  scalar D10;
59  scalar D11;
60  scalar D12;
61  scalar D20;
62  scalar D21;
63  scalar D22;
64  scalar w1;
65  scalar w2;
66 
67  // SA-fv3 coefficients
68  scalar sigmaNut;
69  scalar kappa;
70  scalar Cb1;
71  scalar Cb2;
72  scalar Cw1;
73  scalar Cw2;
74  scalar Cw3;
75  scalar Cv1;
76  scalar Cv2;
77 
78  // SA-fv3 functions
79  tmp<volScalarField> chi(
80  volScalarField& nuTilda,
81  const volScalarField& nu);
82 
83  tmp<volScalarField> fv1(
84  volScalarField& chi);
85 
86  tmp<volScalarField> fv2(
87  volScalarField& chi,
88  volScalarField& fv1);
89 
90  tmp<volScalarField> fv3(
91  volScalarField& chi,
92  volScalarField& fv1);
93 
94  tmp<volScalarField> fw(
95  volScalarField& Stilda,
96  volScalarField& nuTilda,
97  volScalarField& y);
98 
99  tmp<volScalarField> DnuTildaEff(
100  volScalarField& nuTilda,
101  const volScalarField& nu);
102 
103  void correctNut(
104  volScalarField& nut,
105  volScalarField& nuTilda,
106  const volScalarField& nu);
107 
108  // Preimal residuals: meanFlow and SA
109  void calcPriResIrkOrig(
110  volVectorField& U0, // oldTime U
111  volVectorField& U1, // 1st stage
112  volScalarField& p1,
113  surfaceScalarField& phi1,
114  volScalarField& nuTilda1,
115  volScalarField& nut1,
116  volVectorField& U2, // 2nd stage
117  volScalarField& p2,
118  surfaceScalarField& phi2,
119  volScalarField& nuTilda2,
120  volScalarField& nut2,
121  const volScalarField& nu,
122  const scalar& deltaT, // current dt
123  volVectorField& U1Res, // Residual for 1st stage
124  volScalarField& p1Res,
125  surfaceScalarField& phi1Res,
126  volVectorField& U2Res, // Residual for end stage
127  volScalarField& p2Res,
128  surfaceScalarField& phi2Res,
129  const scalar& relaxUEqn);
130 
131  void calcPriSAResIrkOrig(
132  volScalarField& nuTilda0, // oldTime nuTilda
133  volVectorField& U1, // 1st stage
134  surfaceScalarField& phi1,
135  volScalarField& nuTilda1,
136  volVectorField& U2, // 2nd stage
137  surfaceScalarField& phi2,
138  volScalarField& nuTilda2,
139  volScalarField& y,
140  const volScalarField& nu,
141  const scalar& deltaT, // current dt
142  volScalarField& nuTilda1Res, // Residual for 1st stage
143  volScalarField& nuTilda2Res); // Residual for 2nd stage
144 
145  // Some utilities, move to DAUtility later
146 
147  // Swapping two arrays, used for pseudoEqns
148  void swap(List<scalar>& a, List<scalar>& b);
149 
150  // L2 norm
151  scalar L2norm(const List<scalar>& v);
152  scalar L2norm(const List<scalar>& v, const List<scalar>& V);
153  vector L2norm(const List<vector>& U);
154  vector L2norm(const List<vector>& U, const List<scalar>& V);
155  scalar L2norm(const surfaceScalarField& Phi);
156  scalar L2norm(const surfaceScalarField& Phi, const surfaceScalarField& SfArea);
157 
158  // LInf norm for a vector
159  vector getMaxAbs(const List<vector>& U);
160 
161 public:
162  TypeName("DAIrkPimpleFoam");
163  // Constructors
164 
165  //- Construct from components
167  char* argsAll,
168  PyObject* pyOptions);
169 
170  //- Destructor
172  {
173  }
174 
176  virtual void initSolver();
177 
179  virtual label solvePrimal();
180 };
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace Foam
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 #endif
189 
190 // ************************************************************************* //
Foam::DAIrkPimpleFoam::w1
scalar w1
Definition: DAIrkPimpleFoam.H:64
Foam::DAIrkPimpleFoam::Cb1
scalar Cb1
Definition: DAIrkPimpleFoam.H:70
Foam::DAIrkPimpleFoam::calcPriResIrkOrig
void calcPriResIrkOrig(volVectorField &U0, volVectorField &U1, volScalarField &p1, surfaceScalarField &phi1, volScalarField &nuTilda1, volScalarField &nut1, volVectorField &U2, volScalarField &p2, surfaceScalarField &phi2, volScalarField &nuTilda2, volScalarField &nut2, const volScalarField &nu, const scalar &deltaT, volVectorField &U1Res, volScalarField &p1Res, surfaceScalarField &phi1Res, volVectorField &U2Res, volScalarField &p2Res, surfaceScalarField &phi2Res, const scalar &relaxUEqn)
Definition: DAIrkPimpleFoam.C:75
U
volVectorField & U
Definition: createRefsPimpleDyM.H:7
phi2
phi2
Definition: p2EqnIrkPimple.H:19
Foam::DAIrkPimpleFoam::fw
tmp< volScalarField > fw(volScalarField &Stilda, volScalarField &nuTilda, volScalarField &y)
Definition: mySAModel.H:41
Foam::DAIrkPimpleFoam::D12
scalar D12
Definition: DAIrkPimpleFoam.H:60
Foam::DAIrkPimpleFoam::fv2
tmp< volScalarField > fv2(volScalarField &chi, volScalarField &fv1)
Definition: mySAModel.H:21
Foam::DAIrkPimpleFoam::solvePrimal
virtual label solvePrimal()
solve the primal equations
Definition: DAIrkPimpleFoam.C:318
U2
U2
Definition: p2EqnIrkPimple.H:26
Foam::DAIrkPimpleFoam::L2norm
scalar L2norm(const List< scalar > &v)
Definition: myUtilities.H:11
Foam::DAIrkPimpleFoam::w2
scalar w2
Definition: DAIrkPimpleFoam.H:65
Foam::DASolver
Definition: DASolver.H:54
Foam::DAIrkPimpleFoam::Cb2
scalar Cb2
Definition: DAIrkPimpleFoam.H:71
Foam::DAIrkPimpleFoam::chi
tmp< volScalarField > chi(volScalarField &nuTilda, const volScalarField &nu)
Definition: mySAModel.H:7
Foam::DAIrkPimpleFoam::D10
scalar D10
Definition: DAIrkPimpleFoam.H:58
Foam::DAIrkPimpleFoam::Cw3
scalar Cw3
Definition: DAIrkPimpleFoam.H:74
Foam::DAIrkPimpleFoam::fv1
tmp< volScalarField > fv1(volScalarField &chi)
Definition: mySAModel.H:14
Foam::DAIrkPimpleFoam::correctNut
void correctNut(volScalarField &nut, volScalarField &nuTilda, const volScalarField &nu)
Definition: mySAModel.H:70
Foam::DAIrkPimpleFoam::~DAIrkPimpleFoam
virtual ~DAIrkPimpleFoam()
Definition: DAIrkPimpleFoam.H:171
Foam::DAIrkPimpleFoam::DnuTildaEff
tmp< volScalarField > DnuTildaEff(volScalarField &nuTilda, const volScalarField &nu)
Definition: mySAModel.H:62
Foam::DAIrkPimpleFoam::D22
scalar D22
Definition: DAIrkPimpleFoam.H:63
Foam::DAIrkPimpleFoam::Cv1
scalar Cv1
Definition: DAIrkPimpleFoam.H:75
nuTilda2
nuTilda2
Definition: nuTilda2EqnIrkPimple.H:38
phi1
phi1
Definition: p1EqnIrkPimple.H:19
Foam::DAIrkPimpleFoam::Cw1
scalar Cw1
Definition: DAIrkPimpleFoam.H:72
Foam::DAIrkPimpleFoam::swap
void swap(List< scalar > &a, List< scalar > &b)
Definition: myUtilities.H:3
Foam::DAIrkPimpleFoam::D11
scalar D11
Definition: DAIrkPimpleFoam.H:59
Foam
Definition: checkGeometry.C:32
Foam::DAIrkPimpleFoam::kappa
scalar kappa
Definition: DAIrkPimpleFoam.H:69
Foam::DAIrkPimpleFoam::initSolver
virtual void initSolver()
initialize fields and variables
Definition: DAIrkPimpleFoam.C:309
Foam::DAIrkPimpleFoam::sigmaNut
scalar sigmaNut
Definition: DAIrkPimpleFoam.H:68
Foam::DAIrkPimpleFoam::calcPriSAResIrkOrig
void calcPriSAResIrkOrig(volScalarField &nuTilda0, volVectorField &U1, surfaceScalarField &phi1, volScalarField &nuTilda1, volVectorField &U2, surfaceScalarField &phi2, volScalarField &nuTilda2, volScalarField &y, const volScalarField &nu, const scalar &deltaT, volScalarField &nuTilda1Res, volScalarField &nuTilda2Res)
Definition: DAIrkPimpleFoam.C:210
nuTilda1
nuTilda1
Definition: nuTilda1EqnIrkPimple.H:38
Foam::DAIrkPimpleFoam::D20
scalar D20
Definition: DAIrkPimpleFoam.H:61
p1
p1
Definition: p1EqnIrkPimple.H:24
Foam::DAIrkPimpleFoam::DAIrkPimpleFoam
DAIrkPimpleFoam(char *argsAll, PyObject *pyOptions)
Definition: DAIrkPimpleFoam.C:41
Foam::DAIrkPimpleFoam::Cv2
scalar Cv2
Definition: DAIrkPimpleFoam.H:76
Foam::DAIrkPimpleFoam::D21
scalar D21
Definition: DAIrkPimpleFoam.H:62
Foam::DAIrkPimpleFoam::Cw2
scalar Cw2
Definition: DAIrkPimpleFoam.H:73
U1
U1
Definition: p1EqnIrkPimple.H:26
DASolver.H
Foam::DAIrkPimpleFoam::TypeName
TypeName("DAIrkPimpleFoam")
Foam::DAIrkPimpleFoam::getMaxAbs
vector getMaxAbs(const List< vector > &U)
Definition: myUtilities.H:113
Foam::DAIrkPimpleFoam
Definition: DAIrkPimpleFoam.H:52
Foam::DAIrkPimpleFoam::fv3
tmp< volScalarField > fv3(volScalarField &chi, volScalarField &fv1)
Definition: mySAModel.H:28
p2
p2
Definition: p2EqnIrkPimple.H:24