VoFSolver.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 "VoFSolver.H"
27 #include "localEulerDdtScheme.H"
28 #include "linear.H"
29 #include "fvcDiv.H"
30 #include "fvcMeshPhi.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace solvers
37 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 {
48 }
49 
50 
51 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
52 
54 {
55  if (rAU.valid())
56  {
57  rAU() = 1.0/UEqn.A();
58  }
59  else
60  {
61  rAU = (1.0/UEqn.A()).ptr();
62  }
63 }
64 
65 
67 {
68  if (!(correctPhi || mesh.topoChanging()))
69  {
70  rAU.clear();
71  }
72 }
73 
74 
76 {
77  fluidSolver::correctCoNum(phi);
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
84 (
85  fvMesh& mesh,
86  autoPtr<VoFMixture> mixturePtr
87 )
88 :
89  fluidSolver(mesh),
90 
91  mixturePtr_(mixturePtr),
92  mixture_(mixturePtr_()),
93 
94  divAlphaName("div(phi,alpha)"),
95 
96  U_
97  (
98  IOobject
99  (
100  "U",
101  runTime.name(),
102  mesh,
103  IOobject::MUST_READ,
104  IOobject::AUTO_WRITE
105  ),
106  mesh
107  ),
108 
109  phi_
110  (
111  IOobject
112  (
113  "phi",
114  runTime.name(),
115  mesh,
116  IOobject::READ_IF_PRESENT,
117  IOobject::AUTO_WRITE
118  ),
119  linearInterpolate(U_) & mesh.Sf()
120  ),
121 
122  buoyancy(mesh),
123 
124  p_rgh(buoyancy.p_rgh),
125 
126  rho(mixture_.rho()),
127 
128  rhoPhi
129  (
130  IOobject
131  (
132  "rhoPhi",
133  runTime.name(),
134  mesh,
135  IOobject::NO_READ,
136  IOobject::NO_WRITE
137  ),
138  fvc::interpolate(rho)*phi_
139  ),
140 
141  MRF(mesh),
142 
143  mixture(mixture_),
144  U(U_),
145  phi(phi_)
146 {
148 
149  if (LTS)
150  {
151  Info<< "Using LTS" << endl;
152 
154  (
155  new volScalarField
156  (
157  IOobject
158  (
160  runTime.name(),
161  mesh,
164  ),
165  mesh,
167  extrapolatedCalculatedFvPatchScalarField::typeName
168  )
169  );
170  }
171 }
172 
173 
174 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
175 
177 {}
178 
179 
180 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
181 
183 {
184  const scalar maxAlphaCo =
185  runTime.controlDict().lookup<scalar>("maxAlphaCo");
186 
187  scalar deltaT = fluidSolver::maxDeltaT();
188 
189  if (alphaCoNum > small)
190  {
191  deltaT = min(deltaT, maxAlphaCo/alphaCoNum*runTime.deltaTValue());
192  }
193 
194  return deltaT;
195 }
196 
197 
199 {
200  // Read the controls
201  readControls();
202 
203  if ((mesh.dynamic() || MRF.size()) && !Uf.valid())
204  {
205  Info<< "Constructing face momentum Uf" << endl;
206 
207  // Ensure the U BCs are up-to-date before constructing Uf
208  U_.correctBoundaryConditions();
209 
210  Uf = new surfaceVectorField
211  (
212  IOobject
213  (
214  "Uf",
215  runTime.name(),
216  mesh,
219  ),
220  fvc::interpolate(U_)
221  );
222  }
223 
224  if (transient())
225  {
226  correctCoNum();
227  }
228  else if (LTS)
229  {
230  setRDeltaT();
231  }
232 
234 
235  // Store divU from the previous mesh so that it can be mapped
236  // and used in correctPhi to ensure the corrected phi has the
237  // same divergence
238  if (correctPhi && divergent())
239  {
240  // Construct and register divU for mapping
241  divU = new volScalarField
242  (
243  "divU0",
244  fvc::div(fvc::absolute(phi, U))
245  );
246  }
247 
248  // Update the mesh for topology change, mesh to mesh mapping
249  mesh_.update();
250 }
251 
252 
254 {}
255 
256 
258 {
259  divU.clear();
260 }
261 
262 
263 // ************************************************************************* //
scalar alphaCoNum
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.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
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.
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 base-class for the solution of immiscible fluids using a VOF (volume of fluid) pha...
Definition: VoFSolver.H:73
volScalarField & p_rgh
Reference to the buoyant pressure for buoyant cases.
Definition: VoFSolver.H:107
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
Definition: VoFSolver.C:257
tmp< volScalarField > trDeltaT
Optional LTS reciprocal time-step field.
Definition: VoFSolver.H:141
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
Definition: VoFSolver.C:182
VoFSolver(fvMesh &mesh, autoPtr< VoFMixture >)
Construct from region mesh.
Definition: VoFSolver.C:84
void clearrAU()
Clear the cached rAU is no longer needed.
Definition: VoFSolver.C:66
virtual ~VoFSolver()
Destructor.
Definition: VoFSolver.C:176
virtual void prePredictor()=0
Called at the start of the PIMPLE loop.
Definition: VoFSolver.C:253
void continuityErrors()
Calculate and print the continuity errors.
Definition: VoFSolver.C:45
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Definition: VoFSolver.C:198
virtual void correctCoNum()=0
Correct the cached Courant numbers.
Definition: VoFSolver.C:75
void setrAU(const fvVectorMatrix &UEqn)
Set or update the cached rAU.
Definition: VoFSolver.C:53
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 continuityErrors(const surfaceScalarField &phi)
Calculate and print the continuity errors.
Definition: fluidSolver.C:131
virtual scalar maxDeltaT() const
Return the current maximum time-step for stable solution.
Definition: fluidSolver.C:211
A class for managing temporary objects.
Definition: tmp.H:55
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
IOMRFZoneList MRF(mesh)
Calculate the divergence of the given field.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
volScalarField rAU(1.0/UEqn.A())
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.
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
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
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
SurfaceField< vector > surfaceVectorField
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
correctPhi