twoPhaseSolver.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) 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 "twoPhaseSolver.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace solvers
33 {
35 }
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  fvMesh& mesh,
45 )
46 :
47  VoFSolver(mesh, autoPtr<VoFMixture>(mixturePtr.ptr())),
48 
50 
51  alpha1(mixture.alpha1()),
52  alpha2(mixture.alpha2()),
53 
54  alphaRestart
55  (
57  (
58  IOobject::groupName("alphaPhi", alpha1.group()),
59  runTime.name(),
60  mesh,
61  IOobject::READ_IF_PRESENT,
62  IOobject::AUTO_WRITE
63  ).headerOk()
64  ),
65 
66  alphaPhi1
67  (
68  IOobject
69  (
70  IOobject::groupName("alphaPhi", alpha1.group()),
71  runTime.name(),
72  mesh,
73  IOobject::READ_IF_PRESENT,
74  IOobject::AUTO_WRITE
75  ),
76  phi*fvc::interpolate(alpha1)
77  )
78 {
80 
81  if (alphaRestart)
82  {
83  Info << "Restarting alpha" << endl;
84  }
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
91 {}
92 
93 
94 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
95 
97 {
99 
100  // Do not apply previous time-step mesh compression flux
101  // if the mesh topology changed
102  if (mesh().topoChanged())
103  {
104  talphaPhi1Corr0.clear();
105  }
106 }
107 
108 
110 {
112  alphaPredictor();
113  mixture.correct();
114 }
115 
116 
117 // ************************************************************************* //
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
Two-phase VoF mixture.
Definition: VoFMixture.H:52
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:503
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
virtual void prePredictor()=0
Called at the start of the PIMPLE loop.
Definition: VoFSolver.C:253
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Definition: VoFSolver.C:198
Solver module base-class for for 2 immiscible fluids, with optional mesh motion and mesh topology cha...
virtual void prePredictor()
Called at the start of the PIMPLE loop.
twoPhaseSolver(fvMesh &mesh, autoPtr< twoPhaseVoFMixture >)
Construct from region mesh.
bool alphaRestart
Switch indicating if this is a restart.
volScalarField & alpha1
Reference to the phase1-fraction.
virtual ~twoPhaseSolver()
Destructor.
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Class to represent a VoF mixture.
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
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.
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:111
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47