VoFCavitation.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 "VoFCavitation.H"
27 #include "VoFSolver.H"
29 #include "fvcDdt.H"
30 #include "fvcDiv.H"
31 #include "fvmSup.H"
33 
34 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  namespace fv
39  {
41 
43  (
44  fvModel,
47  );
48  }
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const word& sourceName,
57  const word& modelType,
58  const fvMesh& mesh,
59  const dictionary& dict
60 )
61 :
62  fvModel(sourceName, modelType, mesh, dict),
63 
64  mixture_
65  (
66  mesh.lookupObjectRef<incompressibleTwoPhaseVoFMixture>
67  (
68  "phaseProperties"
69  )
70  ),
71 
72  cavitation_(cavitationModel::New(dict, mixture_))
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 {
80  return
81  {
82  mixture_.alpha1().name(),
83  mixture_.alpha2().name(),
84  "U"
85  };
86 }
87 
88 
90 (
91  const volScalarField& alpha,
92  fvMatrix<scalar>& eqn
93 ) const
94 {
95  if (debug)
96  {
97  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
98  }
99 
100  if (&alpha == &mixture_.alpha1() || &alpha == &mixture_.alpha2())
101  {
102  const solvers::VoFSolver& solver =
103  mesh().lookupObject<solvers::VoFSolver>(solver::typeName);
104 
105  const volScalarField& alpha1 = mixture_.alpha1();
106  const volScalarField& alpha2 = mixture_.alpha2();
107 
108  const dimensionedScalar& rho =
109  &alpha == &alpha1 ? mixture_.rho1() : mixture_.rho2();
110 
111  const scalar s = &alpha == &alpha1 ? +1 : -1;
112 
113  // Volume-fraction linearisation
114  if (&alpha == &eqn.psi())
115  {
116  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
117  (
118  cavitation_->mDot12Alpha()
119  );
120  const volScalarField::Internal vDot1Alpha2(mDot12Alpha[0]/rho);
121  const volScalarField::Internal vDot2Alpha1(mDot12Alpha[1]/rho);
122 
123  eqn +=
124  (&alpha == &alpha1 ? vDot1Alpha2 : vDot2Alpha1)
125  - fvm::Sp(vDot1Alpha2 + vDot2Alpha1, eqn.psi());
126  }
127 
128  // Pressure linearisation
129  else if (&eqn.psi() == &solver.p_rgh)
130  {
132  (
133  cavitation_->mDot12P()
134  );
135  const volScalarField::Internal vDot1P(mDot12P[0]/rho);
136  const volScalarField::Internal vDot2P(mDot12P[1]/rho);
137 
139  mesh().lookupObject<volScalarField>("rho");
140  const volScalarField::Internal& gh =
142 
143  eqn +=
144  fvm::Sp(s*(vDot1P - vDot2P), eqn.psi())
145  + s*(vDot1P - vDot2P)*rho*gh
146  - s*(vDot1P - vDot2P)*cavitation_->pSat();
147  }
148 
149  // Explicit non-linearised value. Probably not used.
150  else
151  {
152  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
153  (
154  cavitation_->mDot12Alpha()
155  );
156  const volScalarField::Internal vDot1(mDot12Alpha[0]*alpha2/rho);
157  const volScalarField::Internal vDot2(mDot12Alpha[1]*alpha1/rho);
158 
159  eqn += s*(vDot1 - vDot2);
160  }
161  }
162  else
163  {
165  << "Support for field " << alpha.name() << " is not implemented"
166  << exit(FatalError);
167  }
168 }
169 
170 
172 (
173  const volScalarField& rho,
174  const volVectorField& U,
175  fvMatrix<vector>& eqn
176 ) const
177 {
178  if (debug)
179  {
180  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
181  }
182 
183  if (U.name() == "U")
184  {
185  const surfaceScalarField& rhoPhi =
186  mesh().lookupObject<surfaceScalarField>("rhoPhi");
187 
188  if (&U == &eqn.psi())
189  {
190  eqn += fvm::Sp(fvc::ddt(rho) + fvc::div(rhoPhi), eqn.psi());
191  }
192  else
193  {
194  eqn += (fvc::ddt(rho) + fvc::div(rhoPhi))*U;
195  }
196  }
197  else
198  {
200  << "Support for field " << U.name() << " is not implemented"
201  << exit(FatalError);
202  }
203 }
204 
205 
207 {}
208 
209 
211 {}
212 
213 
215 {}
216 
217 
219 {
220  return true;
221 }
222 
223 
225 {
226  cavitation_->correct();
227 }
228 
229 
230 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
Abstract base class for cavitation models.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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
Finite volume model abstract base class.
Definition: fvModel.H:60
Cavitation fvModel.
Definition: VoFCavitation.H:98
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: VoFCavitation.C:78
virtual void correct()
Correct the cavitation model.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void addSup(const volScalarField &alpha, fvMatrix< scalar > &eqn) const
Add a source to the phase continuity equation.
Definition: VoFCavitation.C:90
VoFCavitation(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: VoFCavitation.C:55
Class to represent a mixture of two constant density phases.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
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.
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
Base solver module base-class for the solution of immiscible fluids using a VOF (volume of fluid) pha...
Definition: VoFSolver.H:73
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 first temporal derivative.
Calculate the divergence of the given field.
Calculate the matrix for implicit and explicit sources.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
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
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
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
labelList fv(nPoints)
dictionary dict