VoFClouds.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) 2021-2024 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 "VoFClouds.H"
28 #include "fvmSup.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
39 
41  (
42  fvModel,
43  VoFClouds,
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const word& sourceName,
55  const word& modelType,
56  const fvMesh& mesh,
57  const dictionary& dict
58 )
59 :
60  fvModel(sourceName, modelType, mesh, dict),
61  phaseName_(dict.lookup("phase")),
62  carrierPhaseName_(dict.lookup("carrierPhase")),
63  thermo_
64  (
65  mesh.lookupObject<fluidThermo>
66  (
67  IOobject::groupName(physicalProperties::typeName, phaseName_)
68  )
69  ),
70  carrierThermo_
71  (
72  mesh.lookupObject<fluidThermo>
73  (
74  IOobject::groupName(physicalProperties::typeName, carrierPhaseName_)
75  )
76  ),
77  clouds_
78  (
79  dict.lookupOrDefault<wordList>
80  (
81  "clouds",
82  parcelCloudList::defaultCloudNames
83  ),
84  carrierThermo_.rho(),
85  mesh.lookupObject<volVectorField>("U"),
86  mesh.lookupObject<uniformDimensionedVectorField>("g"),
87  carrierThermo_
88  ),
89  curTimeIndex_(-1)
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  return wordList({thermo_.rho()().name(), thermo_.he().name(), "U"});
98 }
99 
100 
102 {
103  if (curTimeIndex_ == mesh().time().timeIndex())
104  {
105  return;
106  }
107 
108  clouds_.evolve();
109 
110  curTimeIndex_ = mesh().time().timeIndex();
111 }
112 
113 
115 (
116  const volScalarField& alpha,
117  const volScalarField& rho,
118  fvMatrix<scalar>& eqn
119 ) const
120 {
121  if (debug)
122  {
123  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
124  }
125 
126  if (&rho == &thermo_.rho()())
127  {
128  eqn += clouds_.Srho();
129  }
130  else
131  {
133  << "Support for field " << rho.name() << " is not implemented"
134  << exit(FatalError);
135  }
136 }
137 
138 
140 (
141  const volScalarField& alpha,
142  const volScalarField& rho,
143  const volScalarField& he,
144  fvMatrix<scalar>& eqn
145 ) const
146 {
147  if (debug)
148  {
149  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
150  }
151 
152  if (&he == &thermo_.he())
153  {
154  eqn += clouds_.Sh(eqn.psi());
155  }
156  else
157  {
159  << "Support for field " << he.name() << " is not implemented"
160  << exit(FatalError);
161  }
162 }
163 
164 
166 (
167  const volScalarField& rho,
168  const volVectorField& U,
169  fvMatrix<vector>& eqn
170 ) const
171 {
172  if (debug)
173  {
174  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
175  }
176 
177  if (U.name() == "U")
178  {
179  eqn += clouds_.SU(eqn.psi());
180  }
181  else
182  {
184  << "Support for field " << U.name() << " is not implemented"
185  << exit(FatalError);
186  }
187 }
188 
189 
191 {
192  // Store the particle positions
193  clouds_.storeGlobalPositions();
194 }
195 
196 
198 {
199  clouds_.topoChange(map);
200 }
201 
202 
204 {
205  clouds_.mapMesh(map);
206 }
207 
208 
210 {
211  clouds_.distribute(map);
212 }
213 
214 
216 {
217  return true;
218 }
219 
220 
221 // ************************************************************************* //
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
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:56
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
Finite volume model abstract base class.
Definition: fvModel.H:60
Lagrangian clouds model for VoF simulations.
Definition: VoFClouds.H:65
virtual bool movePoints()
Update for mesh motion.
Definition: VoFClouds.C:215
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: VoFClouds.C:95
virtual void correct()
Solve the Lagrangian clouds and update the sources.
Definition: VoFClouds.C:101
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: VoFClouds.C:197
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: VoFClouds.C:209
VoFClouds(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: VoFClouds.C:53
virtual void preUpdateMesh()
Prepare for mesh_.update.
Definition: VoFClouds.C:190
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: VoFClouds.C:203
virtual void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn) const
Add explicit contribution to phase continuity.
Definition: VoFClouds.C:115
List of parcel clouds, with the same interface as an individual parcel cloud. This is the object that...
A base class for physical properties.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the matrix for implicit and explicit sources.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
label timeIndex
Definition: getTimeIndex.H:4
thermo he()
labelList fv(nPoints)
dictionary dict