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 "fvmSup.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
38  namespace compressible
39  {
41 
43  (
44  fvModel,
47  );
48  }
49  }
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const word& sourceName,
58  const word& modelType,
59  const fvMesh& mesh,
60  const dictionary& dict
61 )
62 :
63  fvModel(sourceName, modelType, mesh, dict),
64 
65  mixture_
66  (
67  mesh.lookupObjectRef<compressibleTwoPhaseVoFMixture>
68  (
69  "phaseProperties"
70  )
71  ),
72 
73  cavitation_(Foam::compressible::cavitationModel::New(dict, mixture_))
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  return
82  {
83  mixture_.rho1().name(),
84  mixture_.rho2().name()
85  };
86 }
87 
88 
90 (
91  const volScalarField& alpha,
92  const volScalarField& rho,
93  fvMatrix<scalar>& eqn
94 ) const
95 {
96  if (debug)
97  {
98  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
99  }
100 
101  if (&rho == &mixture_.rho1() || &rho == &mixture_.rho2())
102  {
103  const solvers::VoFSolver& solver =
104  mesh().lookupObject<solvers::VoFSolver>(solver::typeName);
105 
106  const volScalarField& alpha1 = mixture_.alpha1();
107  const volScalarField& alpha2 = mixture_.alpha2();
108 
109  const scalar s = &alpha == &alpha1 ? +1 : -1;
110 
111  // Volume-fraction linearisation
112  if (&alpha == &eqn.psi())
113  {
114  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
115  (
116  cavitation_->mDot12Alpha()
117  );
118  const volScalarField::Internal& mDot1Alpha2 = mDot12Alpha[0]();
119  const volScalarField::Internal& mDot2Alpha1 = mDot12Alpha[1]();
120 
121  eqn +=
122  (&alpha == &alpha1 ? mDot1Alpha2 : mDot2Alpha1)
123  - fvm::Sp(mDot1Alpha2 + mDot2Alpha1, eqn.psi());
124  }
125 
126  // Pressure linearisation
127  else if (&eqn.psi() == &solver.p_rgh)
128  {
130  (
131  cavitation_->mDot12P()
132  );
133  const volScalarField::Internal& mDot1P = mDot12P[0];
134  const volScalarField::Internal& mDot2P = mDot12P[1];
135 
137  mesh().lookupObject<volScalarField>("rho");
138  const volScalarField::Internal& gh =
140 
141  eqn +=
142  fvm::Sp(s*(mDot1P - mDot2P), eqn.psi())
143  + s*(mDot1P - mDot2P)*rho*gh
144  - s*(mDot1P*cavitation_->pSat1() - mDot2P*cavitation_->pSat2());
145  }
146 
147  // Explicit non-linearised value. Used in density predictors and
148  // continuity error terms.
149  else
150  {
151  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
152  (
153  cavitation_->mDot12Alpha()
154  );
155  const volScalarField::Internal mDot1(mDot12Alpha[0]*alpha2);
156  const volScalarField::Internal mDot2(mDot12Alpha[1]*alpha1);
157 
158  eqn += s*(mDot1 - mDot2);
159  }
160  }
161  else
162  {
164  << "Support for field " << alpha.name() << " is not implemented"
165  << exit(FatalError);
166  }
167 }
168 
169 
171 {}
172 
173 
175 {}
176 
177 
179 (
180  const polyDistributionMap&
181 )
182 {}
183 
184 
186 {
187  return true;
188 }
189 
190 
192 {
193  cavitation_->correct();
194 }
195 
196 
197 // ************************************************************************* //
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.
Class to represent a mixture of two rhoFluidThermo-based phases.
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
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:79
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, const volScalarField &rho, 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:56
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 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))
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)
defineTypeNameAndDebug(VoFCavitation, 0)
addToRunTimeSelectionTable(fvModel, VoFCavitation, dictionary)
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