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-2025 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 
53  pimple.dict().lookupOrDefault<bool>("momentumPredictor", false);
54 
55  faceMomentum =
56  pimple.dict().lookupOrDefault<Switch>("faceMomentum", false);
57 
59  pimple.dict().lookupOrDefault<Switch>("dragCorrection", false);
60 
62  pimple.dict().lookupOrDefault<int>("nEnergyCorrectors", 1);
63 
65 
66  return true;
67 }
68 
69 
70 void Foam::solvers::multiphaseEuler::correctCoNum()
71 {
72  scalarField sumPhi
73  (
74  fvc::surfaceSum(mag(phi))().primitiveField()
75  );
76 
77  forAll(movingPhases, movingPhasei)
78  {
79  sumPhi = max
80  (
81  sumPhi,
82  fvc::surfaceSum(mag(movingPhases[movingPhasei].phi()))()
83  .primitiveField()
84  );
85  }
86 
87  CoNum_ = 0.5*gMax(sumPhi/mesh.V().primitiveField())*runTime.deltaTValue();
88 
89  const scalar meanCoNum =
90  0.5
91  *(gSum(sumPhi)/gSum(mesh.V().primitiveField()))
92  *runTime.deltaTValue();
93 
94  Info<< "Courant Number mean: " << meanCoNum
95  << " max: " << CoNum << endl;
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
102 :
103  fluidSolver(mesh),
104 
105  predictMomentum
106  (
107  pimple.dict().lookupOrDefault<Switch>("momentumPredictor", false)
108  ),
109 
110  faceMomentum
111  (
112  pimple.dict().lookupOrDefault<Switch>("faceMomentum", false)
113  ),
114 
115  dragCorrection
116  (
117  pimple.dict().lookupOrDefault<Switch>("dragCorrection", false)
118  ),
119 
120  nEnergyCorrectors
121  (
122  pimple.dict().lookupOrDefault<int>("nEnergyCorrectors", 1)
123  ),
124 
125  trDeltaT
126  (
127  LTS
128  ? new volScalarField
129  (
130  IOobject
131  (
132  fv::localEulerDdt::rDeltaTName,
133  runTime.name(),
134  mesh,
135  IOobject::READ_IF_PRESENT,
136  IOobject::AUTO_WRITE
137  ),
138  mesh,
140  extrapolatedCalculatedFvPatchScalarField::typeName
141  )
142  : nullptr
143  ),
144 
145  trDeltaTf
146  (
147  LTS && faceMomentum
148  ? new surfaceScalarField
149  (
150  IOobject
151  (
152  fv::localEulerDdt::rDeltaTfName,
153  runTime.name(),
154  mesh,
155  IOobject::READ_IF_PRESENT,
156  IOobject::AUTO_WRITE
157  ),
158  mesh,
160  )
161  : nullptr
162  ),
163 
164  buoyancy(mesh),
165 
166  fluid_(mesh),
167 
168  phases_(fluid_.phases()),
169 
170  movingPhases_(fluid_.movingPhases()),
171 
172  phi_(fluid_.phi()),
173 
174  momentumTransferSystem_(fluid_),
175 
176  heatTransferSystem_(fluid_),
177 
178  populationBalanceSystem_(fluid_),
179 
180  p_(movingPhases_[0].fluidThermo().p()),
181 
182  p_rgh_(buoyancy.p_rgh),
183 
185  (
186  p_,
187  p_rgh_,
188  pimple.dict(),
189  fluid_.incompressible()
190  ),
191 
192  MRF(fluid_.MRF()),
193 
194  fluid(fluid_),
195  phases(phases_),
196  movingPhases(movingPhases_),
197  momentumTransfer(momentumTransferSystem_),
198  heatTransfer(heatTransferSystem_),
199  p(p_),
200  p_rgh(p_rgh_),
201  phi(phi_)
202 {
203  // Read the controls
204  read();
205 
207 
208  if (transient())
209  {
210  correctCoNum();
211  }
212 }
213 
214 
215 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
216 
218 {}
219 
220 
221 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
222 
224 {
225  if (transient())
226  {
227  correctCoNum();
228  }
229  else if (LTS)
230  {
231  setRDeltaT();
232  }
233 
234  // Store divU from the previous mesh so that it can be
235  // mapped and used in correctPhi to ensure the corrected phi
236  // has the same divergence
237  if (correctPhi || mesh.topoChanging())
238  {
239  // Construct and register divU for mapping
240  divU = new volScalarField
241  (
242  "divU0",
243  fvc::div(fvc::absolute(phi, movingPhases[0].U()))
244  );
245  }
246 
248 
249  // Update the mesh for topology change, mesh to mesh mapping
250  mesh_.update();
251 }
252 
253 
255 {
256  if (pimple.thermophysics() || pimple.flow())
257  {
258  alphaControls.correct(CoNum);
259 
260  fluid_.solve(alphaControls, rAs, momentumTransferSystem_);
261  populationBalanceSystem_.solve();
262 
263  fluid_.correct();
264  populationBalanceSystem_.correct();
265 
266  fluid_.correctContinuityError(populationBalanceSystem_.dmdts());
267  }
268 }
269 
270 
272 {
273  fluid_.predictMomentumTransport();
274 }
275 
276 
278 {
279  // Moved inside the nEnergyCorrectors loop in thermophysicalPredictor()
280  // fluid_.predictThermophysicalTransport();
281 }
282 
283 
285 {
286  fluid_.correctMomentumTransport();
287 }
288 
289 
291 {
292  fluid_.correctThermophysicalTransport();
293 }
294 
295 
297 {
298  divU.clear();
299 }
300 
301 
302 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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:307
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 &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1810
bool topoChanging() const
Does the mesh topology change?
Definition: fvMesh.C:705
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:260
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:447
Provides controls for the pressure reference in closed-volume simulations.
virtual const dictionary & dict() const
Return the solution dictionary.
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:348
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:107
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
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
virtual bool read()
Read controls.
Definition: fluidSolver.C:51
Solver module for steady or transient turbulent flow of compressible fluids with heat-transfer for HV...
Definition: fluid.H:70
Solver module for a system of any number of compressible fluid phases with a common pressure,...
virtual void momentumTransportCorrector()
Correct the momentum transport.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
const volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
virtual void momentumTransportPredictor()
Predict the momentum transport.
int nEnergyCorrectors
Number of energy correctors.
Switch predictMomentum
Momentum equation predictor switch.
virtual void thermophysicalTransportCorrector()
Correct the thermophysical transport.
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.
phaseSystem::alphaControl alphaControls
virtual ~multiphaseEuler()
Destructor.
virtual void thermophysicalTransportPredictor()
Predict thermophysical transport.
multiphaseEuler(fvMesh &mesh)
Construct from region mesh.
virtual bool read()
Read controls.
scalar CoNum
scalar meanCoNum
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
IOMRFZoneList MRF(mesh)
pimpleControl pimple(mesh)
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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.
U
Definition: pEqn.H:72
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< VolInternalField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
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:258
const dimensionSet dimless
messageStream Info
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Type gMax(const FieldField< Field, Type > &f)
correctPhi
labelList fv(nPoints)
dictionary dict
volScalarField & p
void read(const dictionary &dict)
Read the alpha and MULES controls from dict.
Definition: phaseSystem.C:415
Foam::surfaceFields.