isothermalFluid.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) 2022-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 "isothermalFluid.H"
27 #include "localEulerDdtScheme.H"
29 #include "fvcMeshPhi.H"
30 #include "fvcVolumeIntegrate.H"
31 #include "fvcReconstruct.H"
32 #include "fvcDiv.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::isothermalFluid::correctCoNum()
51 {
52  fluidSolver::correctCoNum(rho, phi);
53 }
54 
55 
56 void Foam::solvers::isothermalFluid::continuityErrors()
57 {
59 }
60 
61 
62 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
63 
66 (
68 ) const
69 {
70  if (mesh.moving())
71  {
72  return
73  work
74  + fvc::div
75  (
77  p/rho,
78  "div(phi,(p|rho))"
79  )();
80  }
81  else
82  {
83  return move(work);
84  }
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
91 (
92  fvMesh& mesh,
93  autoPtr<fluidThermo> thermoPtr
94 )
95 :
96  fluidSolver(mesh),
97 
98  thermoPtr_(thermoPtr),
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  dpdt
117  (
118  IOobject
119  (
120  "dpdt",
121  runTime.name(),
122  mesh
123  ),
124  mesh,
126  ),
127 
128  buoyancy(buoyancy::New(mesh)),
129 
130  p_rgh(buoyancy.valid() ? buoyancy->p_rgh : p_),
131 
133  (
134  p_,
135  p_rgh,
136  pimple.dict(),
137  thermo_.incompressible()
138  ),
139 
140  U_
141  (
142  IOobject
143  (
144  "U",
145  runTime.name(),
146  mesh,
147  IOobject::MUST_READ,
148  IOobject::AUTO_WRITE
149  ),
150  mesh
151  ),
152 
153  phi_
154  (
155  IOobject
156  (
157  "phi",
158  runTime.name(),
159  mesh,
160  IOobject::READ_IF_PRESENT,
161  IOobject::AUTO_WRITE
162  ),
163  linearInterpolate(rho_*U_) & mesh.Sf()
164  ),
165 
166  K("K", 0.5*magSqr(U_)),
167 
168  momentumTransport
169  (
170  compressible::momentumTransportModel::New
171  (
172  rho_,
173  U_,
174  phi_,
175  thermo_
176  )
177  ),
178 
179  initialMass(fvc::domainIntegrate(rho_)),
180 
181  MRF(mesh),
182 
183  thermo(thermo_),
184  p(p_),
185  rho(rho_),
186  U(U_),
187  phi(phi_)
188 {
190  momentumTransport->validate();
191 
192  if (buoyancy.valid())
193  {
195  (
196  p_rgh,
197  p_,
198  rho_,
199  U,
200  buoyancy->gh,
201  buoyancy->ghf,
202  buoyancy->pRef,
203  thermo_,
204  pimple.dict()
205  );
206 
208  (
209  IOobject
210  (
211  "netForce",
212  runTime.name(),
213  mesh
214  ),
216  (
218  *mesh.magSf()
219  )
220  );
221  }
222 
223  if (transient())
224  {
225  correctCoNum();
226  }
227  else if (LTS)
228  {
229  Info<< "Using LTS" << endl;
230 
232  (
233  new volScalarField
234  (
235  IOobject
236  (
238  runTime.name(),
239  mesh,
242  ),
243  mesh,
245  extrapolatedCalculatedFvPatchScalarField::typeName
246  )
247  );
248  }
249 }
250 
251 
253 :
254  isothermalFluid(mesh, fluidThermo::New(mesh))
255 {}
256 
257 
258 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
259 
261 {}
262 
263 
264 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
265 
267 {
268  if ((mesh.dynamic() || MRF.size()) && !rhoUf.valid())
269  {
270  Info<< "Constructing face momentum rhoUf" << endl;
271 
272  // Ensure the U BCs are up-to-date before constructing Uf
273  U_.correctBoundaryConditions();
274 
275  rhoUf = new surfaceVectorField
276  (
277  IOobject
278  (
279  "rhoUf",
280  runTime.name(),
281  mesh,
284  ),
286  );
287  }
288 
289  if (transient())
290  {
291  correctCoNum();
292  }
293  else if (LTS)
294  {
295  setRDeltaT();
296  }
297 
298  // Store divrhoU from the previous mesh so that it can be mapped
299  // and used in correctPhi to ensure the corrected phi has the
300  // same divergence
301  if (correctPhi || mesh.topoChanging())
302  {
303  divrhoU = new volScalarField
304  (
305  "divrhoU",
306  fvc::div(fvc::absolute(phi, rho, U))
307  );
308  }
309 
311 
312  // Store momentum to set rhoUf for introduced faces
313  if (mesh.topoChanging())
314  {
315  rhoU = new volVectorField("rhoU", rho*U);
316 
317  for (label i = 1; i <= rhoUf().nOldTimes(false); ++ i)
318  {
319  rhoU().oldTimeRef(i) == rho.oldTime(i)*U.oldTime(i);
320  }
321  }
322 
323  // Update the mesh for topology change, mesh to mesh mapping
324  mesh_.update();
325 }
326 
327 
329 {
330  thermo_.correct();
331 }
332 
333 
335 {
336  while (pimple.correct())
337  {
338  if (buoyancy.valid())
339  {
340  correctBuoyantPressure();
341  }
342  else
343  {
344  correctPressure();
345  }
346  }
347 
348  tUEqn.clear();
349 }
350 
351 
353 {
354  if (pimple.correctTransport())
355  {
356  momentumTransport->correct();
357  }
358 }
359 
360 
362 {
363  divrhoU.clear();
364 
365  if (!mesh.schemes().steady())
366  {
367  rho_ = thermo.rho();
368 
369  // Correct rhoUf with the updated density if the mesh is moving
370  fvc::correctRhoUf(rhoUf, rho, U, phi, MRF);
371  }
372 }
373 
374 
375 // ************************************************************************* //
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const word & name() const
Return const reference to name.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
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).
Provides controls for the pressure reference in closed-volume simulations.
virtual const dictionary & dict() const
Return the solution dictionary.
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
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:107
const Time & runTime
Time.
Definition: solver.H:104
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
Buoyancy related data for the Foam::solvers::isothermalFluid solver module when solving buoyant cases...
Definition: buoyancy.H:70
volScalarField gh
(g & h) - ghRef
Definition: buoyancy.H:95
uniformDimensionedScalarField pRef
Reference pressure.
Definition: buoyancy.H:89
surfaceScalarField ghf
(g & hf) - ghRef
Definition: buoyancy.H:98
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
void continuityErrors(const surfaceScalarField &phi)
Calculate and print the continuity errors.
Definition: fluidSolver.C:146
Solver module for steady or transient turbulent flow of compressible isothermal fluids with optional ...
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
volScalarField rho_
The continuity density field.
volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
const surfaceScalarField & phi
Mass-flux field.
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
tmp< volVectorField > netForce
Momentum equation net force source term.
tmp< volScalarField > trDeltaT
Optional LTS reciprocal time-step field.
fluidThermo & thermo_
Reference to the fluid thermophysical properties.
virtual ~isothermalFluid()
Destructor.
const volVectorField & U
Velocity field.
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
tmp< volScalarField::Internal > pressureWork(const tmp< volScalarField::Internal > &) const
Adds the mesh-motion work to the pressure work term provided.
volScalarField & p_
Reference to the pressure field.
isothermalFluid(fvMesh &mesh, autoPtr< fluidThermo >)
Construct from region mesh and thermophysical properties.
const volScalarField & rho
Reference to the continuity density field.
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
const volScalarField & p
Reference to the pressure field.
autoPtr< compressible::momentumTransportModel > momentumTransport
Pointer to the momentum transport model.
A class for managing temporary objects.
Definition: tmp.H:55
tmp< fvVectorMatrix > tUEqn(fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevTau(U)==fvModels.source(rho, U))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
IOMRFZoneList MRF(mesh)
pimpleControl pimple(mesh)
Calculate the divergence of the given field.
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.
K
Definition: pEqn.H:75
U
Definition: pEqn.H:72
dimensionedScalar initialMass
Definition: createFields.H:68
bool valid(const PtrList< ModelType > &l)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
void correctRhoUf(autoPtr< surfaceVectorField > &rhoUf, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const MRFType &MRF)
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:34
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
dimensioned< Type > domainIntegrate(const VolField< Type > &vf)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
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
messageStream Info
void hydrostaticInitialisation(volScalarField &p_rgh, volScalarField &p, volScalarField &rho, const volVectorField &U, const volScalarField &gh, const surfaceScalarField &ghf, const uniformDimensionedScalarField &pRef, fluidThermo &thermo, const dictionary &dict)
tmp< SurfaceField< Type > > linearInterpolate(const VolField< Type > &vf)
Definition: linear.H:108
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
SurfaceField< vector > surfaceVectorField
dimensioned< scalar > magSqr(const dimensioned< Type > &)
correctPhi
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31