shockFluid.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "shockFluid.H"
27 #include "localEulerDdtScheme.H"
29 #include "fvcMeshPhi.H"
30 #include "fvcVolumeIntegrate.H"
31 #include "fvcReconstruct.H"
32 #include "fvcSnGrad.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace solvers
40 {
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::solvers::shockFluid::correctCoNum(const surfaceScalarField& amaxSf)
50 {
51  const scalarField sumAmaxSf(fvc::surfaceSum(amaxSf)().primitiveField());
52 
53  CoNum = 0.5*gMax(sumAmaxSf/mesh.V().field())*runTime.deltaTValue();
54 
55  const scalar meanCoNum =
56  0.5*(gSum(sumAmaxSf)/gSum(mesh.V().field()))*runTime.deltaTValue();
57 
58  Info<< "Courant Number mean: " << meanCoNum
59  << " max: " << CoNum << endl;
60 }
61 
62 
63 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
64 
65 void Foam::solvers::shockFluid::clearTemporaryFields()
66 {
67  rho_pos.clear();
68  rho_neg.clear();
69 
70  rhoU_pos.clear();
71  rhoU_neg.clear();
72 
73  U_pos.clear();
74  U_neg.clear();
75 
76  p_pos.clear();
77  p_neg.clear();
78 
79  a_pos.clear();
80  a_neg.clear();
81 
82  aSf.clear();
83 
84  aphiv_pos.clear();
85  aphiv_neg.clear();
86 
87  devTau.clear();
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
94 :
95  fluidSolver(mesh),
96 
97  thermoPtr_(psiThermo::New(mesh)),
98 
99  thermo_(thermoPtr_()),
100 
101  p_(thermo_.p()),
102 
103  rho_
104  (
105  IOobject
106  (
107  "rho",
108  runTime.name(),
109  mesh,
110  IOobject::READ_IF_PRESENT,
111  IOobject::AUTO_WRITE
112  ),
113  thermo_.renameRho()
114  ),
115 
116  U_
117  (
118  IOobject
119  (
120  "U",
121  runTime.name(),
122  mesh,
123  IOobject::MUST_READ,
124  IOobject::AUTO_WRITE
125  ),
126  mesh
127  ),
128 
129  phi_
130  (
131  IOobject
132  (
133  "phi",
134  runTime.name(),
135  mesh,
136  IOobject::READ_IF_PRESENT,
137  IOobject::AUTO_WRITE
138  ),
139  linearInterpolate(rho_*U_) & mesh.Sf()
140  ),
141 
142  K("K", 0.5*magSqr(U_)),
143 
144  inviscid
145  (
146  max(thermo_.mu()().primitiveField()) > 0
147  ? false
148  : true
149  ),
150 
151  momentumTransport
152  (
153  inviscid
155  : compressible::momentumTransportModel::New
156  (
157  rho_,
158  U_,
159  phi_,
160  thermo_
161  )
162  ),
163 
165  (
166  inviscid
169  (
170  momentumTransport(),
171  thermo_
172  )
173  ),
174 
175  fluxScheme
176  (
177  mesh.schemes().dict().lookupOrDefault<word>("fluxScheme", "Kurganov")
178  ),
179 
180  thermo(thermo_),
181  p(p_),
182  rho(rho_),
183  U(U_),
184  phi(phi_)
185 {
186  // Read the controls
187  readControls();
188 
189  thermo.validate(type(), "e");
190 
191  if (momentumTransport.valid())
192  {
193  momentumTransport->validate();
195  }
196 
197  fluxPredictor();
198 
199  if (transient())
200  {
201  const surfaceScalarField amaxSf
202  (
203  max(mag(aphiv_pos()), mag(aphiv_neg()))
204  );
205 
206  correctCoNum(amaxSf);
207  }
208  else if (LTS)
209  {
210  Info<< "Using LTS" << endl;
211 
213  (
214  new volScalarField
215  (
216  IOobject
217  (
219  runTime.name(),
220  mesh,
223  ),
224  mesh,
226  extrapolatedCalculatedFvPatchScalarField::typeName
227  )
228  );
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
234 
236 {}
237 
238 
239 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
240 
242 {
243  // Read the controls
244  readControls();
245 
246  {
247  const surfaceScalarField amaxSf
248  (
249  max(mag(aphiv_pos()), mag(aphiv_neg()))
250  );
251 
252  if (transient())
253  {
254  correctCoNum(amaxSf);
255  }
256  else if (LTS)
257  {
258  setRDeltaT(amaxSf);
259  }
260  }
261 
263 
264  if (mesh.topoChanging())
265  {
266  pos.clear();
267  neg.clear();
268 
269  clearTemporaryFields();
270  }
271 
272  // Update the mesh for topology change, mesh to mesh mapping
273  mesh_.update();
274 }
275 
276 
278 {
279  if (!inviscid && pimple.correctTransport())
280  {
281  momentumTransport->correct();
282  thermophysicalTransport->correct();
283  }
284 }
285 
286 
288 {}
289 
290 
291 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
Templated abstract base class for thermophysical transport models.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:34
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
Definition: basicThermo.C:331
Base class for single-phase compressible turbulence models.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:260
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:503
static word rDeltaTName
Name of the reciprocal local time-step field.
Abstract base class for turbulence models (RAS, LES and laminar).
Base-class for fluid thermodynamic properties based on compressibility.
Definition: psiThermo.H:54
Abstract base class for run-time selectable region solvers.
Definition: solver.H:55
bool LTS
Switch for local time step transient operation.
Definition: solver.H:69
const Time & runTime
Time.
Definition: solver.H:97
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
void readControls()
Read controls.
Definition: fluidSolver.C:45
scalar CoNum
Current maximum Courant number for time-step control.
Definition: fluidSolver.H:98
Solver module for density-based solution of compressible flow.
Definition: shockFluid.H:73
shockFluid(fvMesh &mesh)
Construct from region mesh.
Definition: shockFluid.C:93
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
Definition: shockFluid.C:287
tmp< surfaceScalarField > aphiv_pos
Definition: shockFluid.H:146
virtual ~shockFluid()
Destructor.
Definition: shockFluid.C:235
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
const volVectorField & U
Reference to the velocity field.
Definition: shockFluid.H:215
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
Definition: shockFluid.C:277
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Definition: shockFluid.C:241
tmp< surfaceScalarField > aphiv_neg
Definition: shockFluid.H:147
autoPtr< compressible::momentumTransportModel > momentumTransport
Pointer to the momentum transport model.
Definition: shockFluid.H:110
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
scalar meanCoNum
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Reconstruct volField from a face flux field.
Calculate the snGrad of the given volField.
Volume integrate volField creating a volField.
Info<< "Creating thermophysical transport model\n"<< endl;turbulenceThermophysicalTransportModels::unityLewisEddyDiffusivity< RASThermophysicalTransportModel< ThermophysicalTransportModel< compressibleMomentumTransportModel, fluidThermo > >> thermophysicalTransport(turbulence(), thermo, true)
K
Definition: pEqn.H:75
U
Definition: pEqn.H:72
const dimensionedScalar mu
Atomic mass unit.
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
Type gSum(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
messageStream Info
tmp< SurfaceField< Type > > linearInterpolate(const VolField< Type > &vf)
Definition: linear.H:108
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTime
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31