compressibleVoF.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 "compressibleVoF.H"
27 #include "localEulerDdtScheme.H"
28 #include "fvcDdt.H"
29 #include "fvcDiv.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace solvers
37 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
45 
47 {
49 
50  const dictionary& alphaControls = mesh.solution().solverDict(alpha1.name());
51 
53  alphaControls.lookupOrDefault("vDotResidualAlpha", 1e-4);
54 
55  return true;
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
62 :
64  (
65  mesh,
67  ),
68 
69  mixture_
70  (
72  ),
73 
74  p(mixture_.p()),
75 
76  vDot
77  (
78  IOobject
79  (
80  "vDot",
81  runTime.name(),
82  mesh,
83  IOobject::READ_IF_PRESENT,
84  IOobject::AUTO_WRITE
85  ),
86  alpha1()*fvc::div(phi)()()
87  ),
88 
89  pressureReference_
90  (
91  p,
92  p_rgh,
93  pimple.dict(),
94  false
95  ),
96 
97  alphaRhoPhi1
98  (
99  IOobject::groupName("alphaRhoPhi", alpha1.group()),
100  fvc::interpolate(mixture_.thermo1().rho())*alphaPhi1
101  ),
102 
103  alphaRhoPhi2
104  (
105  IOobject::groupName("alphaRhoPhi", alpha2.group()),
106  fvc::interpolate(mixture_.thermo2().rho())*alphaPhi2
107  ),
108 
109  K("K", 0.5*magSqr(U)),
110 
111  momentumTransport
112  (
113  rho,
114  U,
115  phi,
116  rhoPhi,
117  alphaPhi1,
118  alphaPhi2,
119  alphaRhoPhi1,
120  alphaRhoPhi2,
121  mixture_
122  ),
123 
124  thermophysicalTransport(momentumTransport),
125 
126  mixture(mixture_)
127 {
128  read();
129 
130  if (correctPhi || mesh.topoChanging())
131  {
132  rAU = new volScalarField
133  (
134  IOobject
135  (
136  "rAU",
137  runTime.name(),
138  mesh,
141  ),
142  mesh,
144  );
145  }
146 }
147 
148 
149 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
150 
152 {}
153 
154 
155 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
156 
158 {
160 
161  const volScalarField& rho1 = mixture_.thermo1().rho();
162  const volScalarField& rho2 = mixture_.thermo2().rho();
163 
164  alphaRhoPhi1 = fvc::interpolate(rho1)*alphaPhi1;
165  alphaRhoPhi2 = fvc::interpolate(rho2)*alphaPhi2;
166 
167  rhoPhi = alphaRhoPhi1 + alphaRhoPhi2;
168 
169  contErr1 =
170  (
171  fvc::ddt(alpha1, rho1)()() + fvc::div(alphaRhoPhi1)()()
172  - (fvModels().source(alpha1, rho1)&rho1)()
173  );
174 
175  contErr2 =
176  (
177  fvc::ddt(alpha2, rho2)()() + fvc::div(alphaRhoPhi2)()()
178  - (fvModels().source(alpha2, rho2)&rho2)()
179  );
180 }
181 
182 
184 {
185  momentumTransport.predict();
186 }
187 
188 
190 {
191  thermophysicalTransport.predict();
192 }
193 
194 
196 {
197  momentumTransport.correct();
198 }
199 
200 
202 {
203  thermophysicalTransport.correct();
204 }
205 
206 
207 // ************************************************************************* //
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Class to represent a mixture of two rhoFluidThermo-based phases.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1810
bool topoChanging() const
Does the mesh topology change?
Definition: fvMesh.C:705
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
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
const Time & runTime
Time.
Definition: solver.H:104
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
autoPtr< volScalarField > rAU
Inverse momentum equation diagonal.
Definition: VoFSolver.H:129
Solver module for 2 compressible, non-isothermal immiscible fluids using a VOF (volume of fluid) phas...
virtual void momentumTransportCorrector()
Correct the momentum transport.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
virtual void momentumTransportPredictor()
Predict the momentum transport.
scalar vDotResidualAlpha
Compressibility source stabilisation tolerance.
compressibleVoF(fvMesh &mesh)
Construct from region mesh.
virtual void thermophysicalTransportCorrector()
Correct the thermophysical transport.
virtual ~compressibleVoF()
Destructor.
virtual void thermophysicalTransportPredictor()
Predict thermophysical transport.
virtual bool read()
Read controls.
bool correctPhi
Switch to correct the flux after mesh change.
Definition: fluidSolver.H:101
virtual void prePredictor()
Called at the start of the PIMPLE loop.
volScalarField & alpha1
Reference to the phase1-fraction.
virtual bool read()
Read controls.
Solver module base-class for 2 immiscible fluids using a VOF (volume of fluid) phase-fraction based i...
Class to represent a VoF mixture.
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Calculate the first temporal derivative.
Calculate the divergence of the given field.
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 char *const group
Group name for atomic constants.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:134
const dimensionSet dimTime
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
dictionary dict
volScalarField & p