incompressibleFluid.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-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 "incompressibleFluid.H"
27 #include "localEulerDdtScheme.H"
28 #include "linear.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace solvers
36 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 {
47  fluidSolver::correctCoNum(phi);
48 }
49 
50 
52 {
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
60 :
61  fluidSolver(mesh),
62 
63  p_
64  (
65  IOobject
66  (
67  "p",
68  runTime.name(),
69  mesh,
70  IOobject::MUST_READ,
71  IOobject::AUTO_WRITE
72  ),
73  mesh
74  ),
75 
77 
78  U_
79  (
80  IOobject
81  (
82  "U",
83  runTime.name(),
84  mesh,
85  IOobject::MUST_READ,
86  IOobject::AUTO_WRITE
87  ),
88  mesh
89  ),
90 
91  phi_
92  (
93  IOobject
94  (
95  "phi",
96  runTime.name(),
97  mesh,
98  IOobject::READ_IF_PRESENT,
99  IOobject::AUTO_WRITE
100  ),
101  linearInterpolate(U_) & mesh.Sf()
102  ),
103 
104  viscosity(viscosityModel::New(mesh)),
105 
106  momentumTransport
107  (
108  incompressible::momentumTransportModel::New
109  (
110  U_,
111  phi_,
112  viscosity
113  )
114  ),
115 
116  MRF(mesh),
117 
118  p(p_),
119  U(U_),
120  phi(phi_)
121 {
123 
124  momentumTransport->validate();
125 
126  if (transient())
127  {
128  correctCoNum();
129  }
130  else if (LTS)
131  {
132  Info<< "Using LTS" << endl;
133 
135  (
136  new volScalarField
137  (
138  IOobject
139  (
141  runTime.name(),
142  mesh,
145  ),
146  mesh,
148  extrapolatedCalculatedFvPatchScalarField::typeName
149  )
150  );
151  }
152 }
153 
154 
155 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
156 
158 {}
159 
160 
161 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
162 
164 {
165  // Read the controls
166  readControls();
167 
168  if ((mesh.dynamic() || MRF.size()) && !Uf.valid())
169  {
170  Info<< "Constructing face momentum Uf" << endl;
171 
172  // Ensure the U BCs are up-to-date before constructing Uf
173  U_.correctBoundaryConditions();
174 
175  Uf = new surfaceVectorField
176  (
177  IOobject
178  (
179  "Uf",
180  runTime.name(),
181  mesh,
184  ),
186  );
187  }
188 
190 
191  if (transient())
192  {
193  correctCoNum();
194  }
195  else if (LTS)
196  {
197  setRDeltaT();
198  }
199 
200  // Update the mesh for topology change, mesh to mesh mapping
201  mesh_.update();
202 }
203 
204 
206 {
207  if (pimple.predictTransport())
208  {
209  momentumTransport->predict();
210  }
211 }
212 
213 
215 {}
216 
217 
219 {
220  while (pimple.correct())
221  {
222  correctPressure();
223  }
224 
225  tUEqn.clear();
226 }
227 
228 
230 {
231  if (pimple.correctTransport())
232  {
233  viscosity->correct();
234  momentumTransport->correct();
235  }
236 }
237 
238 
240 {}
241 
242 
243 // ************************************************************************* //
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
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
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).
Provides controls for the pressure reference in closed-volume simulations.
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 continuityErrors(const surfaceScalarField &phi)
Calculate and print the continuity errors.
Definition: fluidSolver.C:131
Solver module for steady or transient turbulent flow of incompressible isothermal fluids with optiona...
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
const surfaceScalarField & phi
Reference to the volumetric-flux field.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
tmp< volScalarField > trDeltaT
Optional LTS reciprocal time-step field.
autoPtr< incompressible::momentumTransportModel > momentumTransport
Pointer to the momentum transport model.
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
void continuityErrors()
Calculate and print the continuity errors.
void correctCoNum()
Correct the cached Courant numbers.
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
const volScalarField & p
Reference to the pressure field.
incompressibleFluid(fvMesh &mesh)
Construct from region mesh.
A class for managing temporary objects.
Definition: tmp.H:55
An abstract base class for Newtonian viscosity models.
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
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)
U
Definition: pEqn.H:72
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
messageStream Info
tmp< SurfaceField< Type > > linearInterpolate(const VolField< Type > &vf)
Definition: linear.H:108
const dimensionSet dimTime
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
SurfaceField< vector > surfaceVectorField
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dictionary dict
volScalarField & p