multiphaseEuler.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 "multiphaseEuler.H"
27 #include "localEulerDdtScheme.H"
28 #include "surfaceFields.H"
29 #include "fvcDiv.H"
30 #include "fvcSurfaceIntegrate.H"
31 #include "fvcMeshPhi.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace solvers
39 {
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
49 {
51 
52  faceMomentum =
53  pimple.dict().lookupOrDefault<Switch>("faceMomentum", false);
54 
56  pimple.dict().lookupOrDefault<Switch>("dragCorrection", false);
57 
59  pimple.dict().lookupOrDefault<Switch>("partialElimination", false);
60 
62  pimple.dict().lookupOrDefault<int>("nEnergyCorrectors", 1);
63 }
64 
65 
66 void Foam::solvers::multiphaseEuler::correctCoNum()
67 {
68  scalarField sumPhi
69  (
70  fvc::surfaceSum(mag(phi))().primitiveField()
71  );
72 
73  forAll(phases, phasei)
74  {
75  sumPhi = max
76  (
77  sumPhi,
78  fvc::surfaceSum(mag(phases[phasei].phi()))().primitiveField()
79  );
80  }
81 
82  CoNum = 0.5*gMax(sumPhi/mesh.V().field())*runTime.deltaTValue();
83 
84  const scalar meanCoNum =
85  0.5*(gSum(sumPhi)/gSum(mesh.V().field()))*runTime.deltaTValue();
86 
87  Info<< "Courant Number mean: " << meanCoNum
88  << " max: " << CoNum << endl;
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93 
95 :
96  fluidSolver(mesh),
97 
98  faceMomentum
99  (
100  pimple.dict().lookupOrDefault<Switch>("faceMomentum", false)
101  ),
102 
103  dragCorrection
104  (
105  pimple.dict().lookupOrDefault<Switch>("dragCorrection", false)
106  ),
107 
108  partialElimination
109  (
110  pimple.dict().lookupOrDefault<Switch>("partialElimination", false)
111  ),
112 
113  nEnergyCorrectors
114  (
115  pimple.dict().lookupOrDefault<int>("nEnergyCorrectors", 1)
116  ),
117 
118  trDeltaT
119  (
120  LTS
121  ? new volScalarField
122  (
123  IOobject
124  (
125  fv::localEulerDdt::rDeltaTName,
126  runTime.name(),
127  mesh,
128  IOobject::READ_IF_PRESENT,
129  IOobject::AUTO_WRITE
130  ),
131  mesh,
133  extrapolatedCalculatedFvPatchScalarField::typeName
134  )
135  : nullptr
136  ),
137 
138  trDeltaTf
139  (
140  LTS && faceMomentum
141  ? new surfaceScalarField
142  (
143  IOobject
144  (
145  fv::localEulerDdt::rDeltaTfName,
146  runTime.name(),
147  mesh,
148  IOobject::READ_IF_PRESENT,
149  IOobject::AUTO_WRITE
150  ),
151  mesh,
153  )
154  : nullptr
155  ),
156 
157  buoyancy(mesh),
158 
159  fluidPtr_(phaseSystem::New(mesh)),
160 
161  fluid_(fluidPtr_()),
162 
163  phases_(fluid_.phases()),
164 
165  phi_(fluid_.phi()),
166 
167  p_(phases_[0].thermo().p()),
168 
169  p_rgh(buoyancy.p_rgh),
170 
172  (
173  p_,
174  p_rgh,
175  pimple.dict(),
176  fluid_.incompressible()
177  ),
178 
179  MRF(fluid_.MRF()),
180 
181  fluid(fluid_),
182  phases(phases_),
183  p(p_),
184  phi(phi_)
185 {
186  // Read the controls
187  readControls();
188 
190 
191  if (transient())
192  {
193  correctCoNum();
194  }
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
199 
201 {}
202 
203 
204 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
205 
207 {
208  // Read the controls
209  readControls();
210 
211  if (transient())
212  {
213  correctCoNum();
214  }
215  else if (LTS)
216  {
217  setRDeltaT();
218  }
219 
220  // Store divU from the previous mesh so that it can be
221  // mapped and used in correctPhi to ensure the corrected phi
222  // has the same divergence
223  if (correctPhi || mesh.topoChanging())
224  {
225  // Construct and register divU for mapping
226  divU = new volScalarField
227  (
228  "divU0",
229  fvc::div
230  (
231  fvc::absolute(phi, fluid.movingPhases()[0].U())
232  )
233  );
234  }
235 
237 
238  // Update the mesh for topology change, mesh to mesh mapping
239  mesh_.update();
240 }
241 
242 
244 {
245  if (pimple.thermophysics() || pimple.flow())
246  {
247  fluid_.solve(rAUs, rAUfs);
248  fluid_.correct();
249  fluid_.correctContinuityError();
250  }
251 
252  if (pimple.flow() && pimple.predictTransport())
253  {
254  fluid_.predictMomentumTransport();
255  }
256 }
257 
258 
260 {
261  if (pimple.flow() && pimple.correctTransport())
262  {
263  fluid_.correctMomentumTransport();
264  fluid_.correctThermophysicalTransport();
265  }
266 }
267 
268 
270 {
271  divU.clear();
272 }
273 
274 
275 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
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:55
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:100
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
Buoyancy related data for the Foam::solvers::isothermalFluid solver module when solving buoyant cases...
Definition: buoyancy.H:70
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
void readControls()
Read controls.
Definition: fluidSolver.C:45
Solver module for steady or transient turbulent flow of compressible fluids with heat-transfer for HV...
Definition: fluid.H:70
const volVectorField & U
Velocity field.
Solver module for a system of any number of compressible fluid phases with a common pressure,...
volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
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.
int nEnergyCorrectors
Number of energy correctors.
Switch partialElimination
Partial elimination drag contribution optimisation.
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
void readControls()
Read controls.
Switch faceMomentum
Cell/face momentum equation switch.
Switch dragCorrection
Cell/face drag correction for cell momentum corrector.
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
virtual ~multiphaseEuler()
Destructor.
multiphaseEuler(fvMesh &mesh)
Construct from region mesh.
scalar CoNum
scalar meanCoNum
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
IOMRFZoneList MRF(mesh)
pimpleControl pimple(mesh)
tmp< surfaceScalarField > trDeltaTf
Definition: createRDeltaTf.H:1
Calculate the divergence of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
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
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
messageStream Info
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
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)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Type gMax(const FieldField< Field, Type > &f)
correctPhi
labelList fv(nPoints)
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31
Foam::surfaceFields.