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-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 "VoFCavitation.H"
28 #include "fvcDdt.H"
29 #include "fvcDiv.H"
30 #include "fvmSup.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  namespace fv
38  {
40 
42  (
43  fvModel,
46  );
47  }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const word& sourceName,
56  const word& modelType,
57  const fvMesh& mesh,
58  const dictionary& dict
59 )
60 :
61  fvModel(sourceName, modelType, mesh, dict),
62 
63  mixture_
64  (
65  mesh.lookupObjectRef<incompressibleTwoPhaseVoFMixture>
66  (
67  "phaseProperties"
68  )
69  ),
70 
71  cavitation_(cavitationModel::New(dict, mixture_))
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
78 {
79  return
80  {
81  mixture_.alpha1().name(),
82  mixture_.alpha2().name(),
83  "U"
84  };
85 }
86 
87 
89 (
90  const volScalarField& alpha,
91  fvMatrix<scalar>& eqn
92 ) const
93 {
94  if (debug)
95  {
96  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
97  }
98 
99  if (&alpha == &mixture_.alpha1() || &alpha == &mixture_.alpha2())
100  {
101  const volScalarField& alpha1 = mixture_.alpha1();
102  const volScalarField& alpha2 = mixture_.alpha2();
103 
104  const dimensionedScalar& rho =
105  &alpha == &alpha1 ? mixture_.rho1() : mixture_.rho2();
106 
107  const scalar s = &alpha == &alpha1 ? +1 : -1;
108 
109  // Volume-fraction linearisation
110  if (&alpha == &eqn.psi())
111  {
112  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
113  (
114  cavitation_->mDot12Alpha()
115  );
116  const volScalarField::Internal vDot1Alpha2(mDot12Alpha[0]/rho);
117  const volScalarField::Internal vDot2Alpha1(mDot12Alpha[1]/rho);
118 
119  eqn +=
120  (&alpha == &alpha1 ? vDot1Alpha2 : vDot2Alpha1)
121  - fvm::Sp(vDot1Alpha2 + vDot2Alpha1, eqn.psi());
122  }
123 
124  // Pressure linearisation
125  else if (eqn.psi().name() == "p_rgh")
126  {
128  (
129  cavitation_->mDot12P()
130  );
131  const volScalarField::Internal vDot1P(mDot12P[0]/rho);
132  const volScalarField::Internal vDot2P(mDot12P[1]/rho);
133 
135  mesh().lookupObject<volScalarField>("rho");
136  const volScalarField::Internal& gh =
137  mesh().lookupObject<volScalarField>("gh");
138 
139  eqn +=
140  fvm::Sp(s*(vDot1P - vDot2P), eqn.psi())
141  + s*(vDot1P - vDot2P)*rho*gh
142  - s*(vDot1P - vDot2P)*cavitation_->pSat();
143  }
144 
145  // Explicit non-linearised value. Probably not used.
146  else
147  {
148  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
149  (
150  cavitation_->mDot12Alpha()
151  );
152  const volScalarField::Internal vDot1(mDot12Alpha[0]*alpha2/rho);
153  const volScalarField::Internal vDot2(mDot12Alpha[1]*alpha1/rho);
154 
155  eqn += s*(vDot1 - vDot2);
156  }
157  }
158  else
159  {
161  << "Support for field " << alpha.name() << " is not implemented"
162  << exit(FatalError);
163  }
164 }
165 
166 
168 (
169  const volScalarField& rho,
170  const volVectorField& U,
171  fvMatrix<vector>& eqn
172 ) const
173 {
174  if (debug)
175  {
176  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
177  }
178 
179  if (U.name() == "U")
180  {
181  const surfaceScalarField& rhoPhi =
182  mesh().lookupObject<surfaceScalarField>("rhoPhi");
183 
184  if (&U == &eqn.psi())
185  {
186  eqn += fvm::Sp(fvc::ddt(rho) + fvc::div(rhoPhi), eqn.psi());
187  }
188  else
189  {
190  eqn += (fvc::ddt(rho) + fvc::div(rhoPhi))*U;
191  }
192  }
193  else
194  {
196  << "Support for field " << U.name() << " is not implemented"
197  << exit(FatalError);
198  }
199 }
200 
201 
203 {}
204 
205 
207 {}
208 
209 
211 {}
212 
213 
215 {
216  return true;
217 }
218 
219 
221 {
222  cavitation_->correct();
223 }
224 
225 
226 // ************************************************************************* //
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 <T> with first() and second() elements.
Definition: Pair.H:65
Abstract base class for cavitation models.
A list of keyword definitions, which are a keyword followed by any number of values (e....
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:99
Finite volume model abstract base class.
Definition: fvModel.H:59
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:77
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:89
VoFCavitation(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: VoFCavitation.C:54
Class to represent a mixture of two constant density phases.
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
#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(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), 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:257
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