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-2024 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 "fvMeshStitcher.H"
28 #include "localEulerDdtScheme.H"
30 #include "fvcMeshPhi.H"
31 #include "fvcVolumeIntegrate.H"
32 #include "fvcReconstruct.H"
33 #include "fvcSnGrad.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace solvers
41 {
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::solvers::shockFluid::correctCoNum(const surfaceScalarField& amaxSf)
51 {
52  const scalarField sumAmaxSf(fvc::surfaceSum(amaxSf)().primitiveField());
53 
54  CoNum_ =
55  0.5*gMax(sumAmaxSf/mesh.V().primitiveField())*runTime.deltaTValue();
56 
57  const scalar meanCoNum =
58  0.5
59  *(gSum(sumAmaxSf)/gSum(mesh.V().primitiveField()))
61 
62  Info<< "Courant Number mean: " << meanCoNum
63  << " max: " << CoNum << endl;
64 }
65 
66 
67 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
68 
69 void Foam::solvers::shockFluid::clearTemporaryFields()
70 {
71  rho_pos.clear();
72  rho_neg.clear();
73 
74  rhoU_pos.clear();
75  rhoU_neg.clear();
76 
77  U_pos.clear();
78  U_neg.clear();
79 
80  p_pos.clear();
81  p_neg.clear();
82 
83  a_pos.clear();
84  a_neg.clear();
85 
86  aSf.clear();
87 
88  aphiv_pos.clear();
89  aphiv_neg.clear();
90 
91  devTau.clear();
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96 
98 :
99  fluidSolver(mesh),
100 
101  thermoPtr_(psiThermo::New(mesh)),
102 
103  thermo_(thermoPtr_()),
104 
105  p_(thermo_.p()),
106 
107  rho_
108  (
109  IOobject
110  (
111  "rho",
112  runTime.name(),
113  mesh,
114  IOobject::READ_IF_PRESENT,
115  IOobject::AUTO_WRITE
116  ),
117  thermo_.renameRho()
118  ),
119 
120  U_
121  (
122  IOobject
123  (
124  "U",
125  runTime.name(),
126  mesh,
127  IOobject::MUST_READ,
128  IOobject::AUTO_WRITE
129  ),
130  mesh
131  ),
132 
133  phi_
134  (
135  IOobject
136  (
137  "phi",
138  runTime.name(),
139  mesh,
140  IOobject::READ_IF_PRESENT,
141  IOobject::AUTO_WRITE
142  ),
143  linearInterpolate(rho_*U_) & mesh.Sf()
144  ),
145 
146  K("K", 0.5*magSqr(U_)),
147 
148  inviscid
149  (
150  max(thermo_.mu().primitiveField()) > 0
151  ? false
152  : true
153  ),
154 
155  momentumTransport
156  (
157  inviscid
159  : compressible::momentumTransportModel::New
160  (
161  rho_,
162  U_,
163  phi_,
164  thermo_
165  )
166  ),
167 
169  (
170  inviscid
173  (
174  momentumTransport(),
175  thermo_
176  )
177  ),
178 
179  fluxScheme
180  (
181  mesh.schemes().dict().lookupOrDefault<word>("fluxScheme", "Kurganov")
182  ),
183 
184  thermo(thermo_),
185  p(p_),
186  rho(rho_),
187  U(U_),
188  phi(phi_)
189 {
190  thermo.validate(type(), "e");
191 
192  if (momentumTransport.valid())
193  {
194  momentumTransport->validate();
196  }
197 
198  fluxPredictor();
199 
200  if (transient())
201  {
202  const surfaceScalarField amaxSf
203  (
204  max(mag(aphiv_pos()), mag(aphiv_neg()))
205  );
206 
207  correctCoNum(amaxSf);
208  }
209  else if (LTS)
210  {
211  Info<< "Using LTS" << endl;
212 
214  (
215  new volScalarField
216  (
217  IOobject
218  (
220  runTime.name(),
221  mesh,
224  ),
225  mesh,
227  extrapolatedCalculatedFvPatchScalarField::typeName
228  )
229  );
230  }
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
235 
237 {}
238 
239 
240 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
241 
243 {
244  {
245  const surfaceScalarField amaxSf
246  (
247  max(mag(aphiv_pos()), mag(aphiv_neg()))
248  );
249 
250  if (transient())
251  {
252  correctCoNum(amaxSf);
253  }
254  else if (LTS)
255  {
256  setRDeltaT(amaxSf);
257  }
258  }
259 
261 
262  if (mesh.topoChanging() || mesh.stitcher().stitches())
263  {
264  pos.clear();
265  neg.clear();
266 
267  clearTemporaryFields();
268  }
269 
270  // Update the mesh for topology change, mesh to mesh mapping
271  mesh_.update();
272 }
273 
274 
276 {
277  if (!inviscid && pimple.correctTransport())
278  {
279  momentumTransport->correct();
280  thermophysicalTransport->correct();
281  }
282 }
283 
284 
286 {}
287 
288 
289 // ************************************************************************* //
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:308
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:99
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
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:56
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
bool LTS
Switch for local time step transient operation.
Definition: solver.H:70
const Time & runTime
Time.
Definition: solver.H:104
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
const scalar & CoNum
Current maximum Courant number for time-step control.
Definition: fluidSolver.H:139
scalar CoNum_
Current maximum Courant number for time-step control.
Definition: fluidSolver.H:104
Solver module for density-based solution of compressible flow.
Definition: shockFluid.H:73
shockFluid(fvMesh &mesh)
Construct from region mesh.
Definition: shockFluid.C:97
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
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
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:275
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
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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)
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