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-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 "twoPhaseSolver.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace solvers
33 {
35 }
36 }
37 
38 
39 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
40 
42 {
44 
45  const dictionary& alphaControls = mesh.solution().solverDict(alpha1.name());
46 
48  (
49  alphaControls.found("nAlphaSubCycles")
50  ? "nAlphaSubCycles"
51  : "nSubCycles",
52  dimless,
53  dimless,
54  alphaControls
55  );
56 
58  (
59  {"nCorrectors", "nAlphaCorr"},
60  1
61  );
62 
63  MULESCorr = alphaControls.lookupOrDefault<Switch>("MULESCorr", false);
64 
66  alphaControls.lookupOrDefault<Switch>("alphaApplyPrevCorr", false);
67 
68  MULEScontrols.read(alphaControls);
69 
70  return true;
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
75 
77 (
78  fvMesh& mesh,
80 )
81 :
82  VoFSolver(mesh, autoPtr<VoFMixture>(mixturePtr.ptr())),
83 
85 
86  alpha1(mixture.alpha1()),
87  alpha2(mixture.alpha2()),
88 
89  alphaRestart
90  (
92  (
93  IOobject::groupName("alphaPhi", alpha1.group()),
94  runTime.name(),
95  mesh,
96  IOobject::READ_IF_PRESENT,
97  IOobject::AUTO_WRITE
98  ).headerOk()
99  ),
100 
101  alphaPhi1
102  (
103  IOobject
104  (
105  IOobject::groupName("alphaPhi", alpha1.group()),
106  runTime.name(),
107  mesh,
108  IOobject::READ_IF_PRESENT,
109  IOobject::AUTO_WRITE
110  ),
111  phi*fvc::interpolate(alpha1)
112  ),
113  alphaPhi2
114  (
115  IOobject
116  (
117  IOobject::groupName("alphaPhi", alpha2.group()),
118  runTime.name(),
119  mesh
120  ),
121  phi - alphaPhi1
122  )
123 {
125 
126  read();
127 
128  if (alphaRestart)
129  {
130  Info << "Restarting alpha" << endl;
131  }
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
136 
138 {}
139 
140 
141 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
142 
144 {
146 
147  // Do not apply previous time-step mesh compression flux
148  // if the mesh topology changed
149  if (mesh().topoChanged())
150  {
151  talphaPhi1Corr0.clear();
152  }
153 }
154 
155 
157 {
158  alphaPredictor();
159  mixture.correct();
160 }
161 
162 
163 // ************************************************************************* //
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
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
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
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
T lookupOrDefaultBackwardsCompatible(const wordList &, const T &) const
Find and return a T, trying a list of keywords in sequence,.
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1799
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1810
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:447
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
Definition: solution.C:348
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
Base solver module base-class for the solution of immiscible fluids using a VOF (volume of fluid) pha...
Definition: VoFSolver.H:73
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
Definition: VoFSolver.C:199
virtual bool read()
Read controls.
Definition: fluidSolver.C:51
Solver module base-class for 2 immiscible fluids, with optional mesh motion and mesh topology changes...
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.
bool alphaApplyPrevCorr
Apply the compression correction from the previous iteration.
MULES::control MULEScontrols
MULES controls.
label nAlphaCorr
Number of alpha correctors.
volScalarField & alpha1
Reference to the phase1-fraction.
autoPtr< Function1< scalar > > nAlphaSubCyclesPtr
Function to calculate the number of explicit MULES sub-cycles.
bool MULESCorr
Semi-implicit MULES switch.
virtual ~twoPhaseSolver()
Destructor.
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
virtual bool read()
Read controls.
Class to represent a VoF mixture.
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:530
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:134
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
messageStream Info
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void read(const dictionary &dict)
Read dict and set the controls.
Definition: MULES.C:44