shockFluid.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::solvers::shockFluid
26 
27 Description
28  Solver module for density-based solution of compressible flow
29 
30  Based on central-upwind schemes of Kurganov and Tadmor with support for
31  mesh-motion and topology change.
32 
33  Reference:
34  \verbatim
35  Greenshields, C. J., Weller, H. G., Gasparini, L.,
36  & Reese, J. M. (2010).
37  Implementation of semi‐discrete, non‐staggered central schemes
38  in a colocated, polyhedral, finite volume framework,
39  for high‐speed viscous flows.
40  International journal for numerical methods in fluids, 63(1), 1-21.
41  \endverbatim
42 
43 SourceFiles
44  shockFluid.C
45 
46 See also
47  Foam::solvers::fluidSolver
48  Foam::solvers::incompressibleFluid
49 
50 \*---------------------------------------------------------------------------*/
51 
52 #ifndef shockFluid_H
53 #define shockFluid_H
54 
55 #include "fluidSolver.H"
56 #include "psiThermo.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace solvers
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class shockFluid Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 class shockFluid
72 :
73  public fluidSolver
74 {
75 
76 protected:
77 
78  // Thermophysical properties
79 
80  //- Pointer to the fluid thermophysical properties
82 
83  //- Reference to the fluid thermophysical properties
85 
86  //- Reference to the pressure field
88 
89  //- The continuity density field
91 
92 
93  // Kinematic properties
94 
95  //- Velocity field
97 
98  //- Mass-flux field
100 
101  //- Kinetic energy field
102  // Used in the energy equation
104 
105 
106  // Momentum transport
107 
108  bool inviscid;
109 
110  //- Pointer to the momentum transport model
112 
113 
114  // Thermophysical transport
115 
118 
119 
120  // Controls
121 
123 
124 
125  // Cached temporary fields
126 
129 
132 
135 
138 
141 
144 
146 
149 
151 
152  //- Optional LTS reciprocal time-step field
154 
155 
156 private:
157 
158  // Private Member Functions
159 
160  //- Set rDeltaT for LTS
161  void setRDeltaT(const surfaceScalarField& amaxSf);
162 
163  //- Correct the cached Courant numbers
164  void correctCoNum(const surfaceScalarField& amaxSf);
165 
166  //- Construct the continuity equation and correct the density
167  void correctDensity();
168 
169  //- Interpolate field vf according to direction dir
170  template<class Type>
171  inline tmp<SurfaceField<Type>> interpolate
172  (
173  const VolField<Type>& vf,
174  const surfaceScalarField& dir,
175  const word& reconFieldName = word::null
176  )
177  {
179  (
181  (
182  vf,
183  dir,
184  "reconstruct("
185  + (reconFieldName != word::null ? reconFieldName : vf.name())
186  + ')'
187  )
188  );
189 
190  SurfaceField<Type>& sf = tsf.ref();
191 
192  sf.rename(vf.name() + '_' + dir.name());
193 
194  return tsf;
195  }
196 
197  void fluxPredictor();
198 
199  void clearTemporaryFields();
200 
201 
202 public:
203 
204  // Public Data
205 
206  //- Reference to the fluid thermophysical properties
207  const psiThermo& thermo;
208 
209  //- Reference to the pressure field
210  const volScalarField& p;
211 
212  //- Reference to the continuity density field
213  const volScalarField& rho;
214 
215  //- Reference to the velocity field
216  const volVectorField& U;
217 
218  //- Reference to the mass-flux field
219  const surfaceScalarField& phi;
220 
221 
222  //- Runtime type information
223  TypeName("shockFluid");
224 
225 
226  // Constructors
227 
228  //- Construct from region mesh
230 
231  //- Disallow default bitwise copy construction
232  shockFluid(const shockFluid&) = delete;
233 
234 
235  //- Destructor
236  virtual ~shockFluid();
237 
238 
239  // Member Functions
240 
241  //- Called at the start of the time-step, before the PIMPLE loop
242  virtual void preSolve();
243 
244  //- Called at the start of the PIMPLE loop to move the mesh
245  virtual void moveMesh();
246 
247  //- Corrections that follow mesh motion
248  virtual void motionCorrector();
249 
250  //- Called at the start of the PIMPLE loop
251  virtual void prePredictor();
252 
253  //- Construct and optionally solve the momentum equation
254  virtual void momentumPredictor();
255 
256  //- Construct and solve the energy equation,
257  // convert to temperature
258  // and update thermophysical and transport properties
259  virtual void thermophysicalPredictor();
260 
261  //- Construct and solve the pressure equation in the PISO loop
262  virtual void pressureCorrector();
263 
264  //- Correct the momentum and thermophysical transport modelling
265  virtual void postCorrector();
266 
267  //- Called after the PIMPLE loop at the end of the time-step
268  virtual void postSolve();
269 
270 
271  // Member Operators
272 
273  //- Disallow default bitwise assignment
274  void operator=(const shockFluid&) = delete;
275 };
276 
277 
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 
280 } // End namespace solvers
281 } // End namespace Foam
282 
283 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
284 
285 #endif
286 
287 // ************************************************************************* //
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Base-class for fluid thermodynamic properties based on compressibility.
Definition: psiThermo.H:56
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
Solver module for density-based solution of compressible flow.
Definition: shockFluid.H:73
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
volScalarField rho_
The continuity density field.
Definition: shockFluid.H:89
tmp< surfaceVectorField > U_neg
Definition: shockFluid.H:136
tmp< surfaceVectorField > U_pos
Definition: shockFluid.H:135
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: shockFluid.H:218
volVectorField U_
Velocity field.
Definition: shockFluid.H:95
tmp< surfaceVectorField > devTau
Definition: shockFluid.H:149
virtual void prePredictor()
Called at the start of the PIMPLE loop.
Definition: prePredictor.C:30
shockFluid(fvMesh &mesh)
Construct from region mesh.
Definition: shockFluid.C:97
tmp< surfaceScalarField > p_neg
Definition: shockFluid.H:139
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
Definition: shockFluid.C:285
tmp< surfaceScalarField > aphiv_pos
Definition: shockFluid.H:146
tmp< surfaceScalarField > neg
Definition: shockFluid.H:127
virtual ~shockFluid()
Destructor.
Definition: shockFluid.C:236
const psiThermo & thermo
Reference to the fluid thermophysical properties.
Definition: shockFluid.H:206
tmp< volScalarField > trDeltaT
Optional LTS reciprocal time-step field.
Definition: shockFluid.H:152
tmp< surfaceVectorField > rhoU_pos
Definition: shockFluid.H:132
volScalarField K
Kinetic energy field.
Definition: shockFluid.H:102
virtual void moveMesh()
Called at the start of the PIMPLE loop to move the mesh.
Definition: moveMesh.C:30
psiThermo & thermo_
Reference to the fluid thermophysical properties.
Definition: shockFluid.H:83
const volVectorField & U
Reference to the velocity field.
Definition: shockFluid.H:215
tmp< surfaceScalarField > p_pos
Definition: shockFluid.H:138
void operator=(const shockFluid &)=delete
Disallow default bitwise assignment.
virtual void motionCorrector()
Corrections that follow mesh motion.
Definition: moveMesh.C:40
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
tmp< surfaceScalarField > rho_neg
Definition: shockFluid.H:130
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
Definition: shockFluid.C:275
tmp< surfaceScalarField > pos
Definition: shockFluid.H:126
tmp< surfaceScalarField > aSf
Definition: shockFluid.H:144
virtual void momentumPredictor()
Construct and optionally solve the momentum equation.
volScalarField & p_
Reference to the pressure field.
Definition: shockFluid.H:86
tmp< surfaceScalarField > a_pos
Definition: shockFluid.H:141
surfaceScalarField phi_
Mass-flux field.
Definition: shockFluid.H:98
autoPtr< fluidThermoThermophysicalTransportModel > thermophysicalTransport
Definition: shockFluid.H:116
const volScalarField & rho
Reference to the continuity density field.
Definition: shockFluid.H:212
tmp< surfaceScalarField > a_neg
Definition: shockFluid.H:142
tmp< surfaceVectorField > rhoU_neg
Definition: shockFluid.H:133
tmp< surfaceScalarField > rho_pos
Definition: shockFluid.H:129
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Definition: shockFluid.C:242
tmp< surfaceScalarField > aphiv_neg
Definition: shockFluid.H:147
const volScalarField & p
Reference to the pressure field.
Definition: shockFluid.H:209
TypeName("shockFluid")
Runtime type information.
autoPtr< compressible::momentumTransportModel > momentumTransport
Pointer to the momentum transport model.
Definition: shockFluid.H:110
autoPtr< psiThermo > thermoPtr_
Pointer to the fluid thermophysical properties.
Definition: shockFluid.H:80
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
volScalarField sf(fieldObject, mesh)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.