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  alphaName_(mixture_.alpha1().name())
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  return {alphaName_, "p_rgh", "U"};
82 }
83 
84 
86 (
87  fvMatrix<scalar>& eqn,
88  const word& fieldName
89 ) const
90 {
91  if (debug)
92  {
93  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
94  }
95 
96  if (fieldName == alphaName_)
97  {
98  const volScalarField::Internal alpha1Coeff
99  (
100  1.0/mixture_.rho1()
101  - mixture_.alpha1()()*(1.0/mixture_.rho1() - 1.0/mixture_.rho2())
102  );
103 
104  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
105  (
106  cavitation_->mDot12Alpha()
107  );
108 
109  const volScalarField::Internal vDot1Alpha(alpha1Coeff*mDot12Alpha[0]());
110  const volScalarField::Internal vDot2Alpha(alpha1Coeff*mDot12Alpha[1]());
111 
112  eqn += fvm::Sp(-vDot2Alpha - vDot1Alpha, eqn.psi()) + vDot1Alpha;
113  }
114 }
115 
116 
118 (
119  const volScalarField&,
120  fvMatrix<scalar>& eqn,
121  const word& fieldName
122 ) const
123 {
124  if (debug)
125  {
126  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
127  }
128 
129  if (fieldName == "p_rgh")
130  {
132  mesh().lookupObject<volScalarField>("rho");
133 
134  const volScalarField::Internal& gh =
135  mesh().lookupObject<volScalarField>("gh");
136 
137  const dimensionedScalar pCoeff
138  (
139  1.0/mixture_.rho1() - 1.0/mixture_.rho2()
140  );
141 
143  (
144  cavitation_->mDot12P()
145  );
146 
147  const volScalarField::Internal vDot1P(pCoeff*mDot12P[0]);
148  const volScalarField::Internal vDot2P(pCoeff*mDot12P[1]);
149 
150  eqn +=
151  (vDot2P - vDot1P)*(cavitation_->pSat() - rho*gh)
152  - fvm::Sp(vDot2P - vDot1P, eqn.psi());
153  }
154 }
155 
156 
158 (
159  const volScalarField& rho,
160  fvMatrix<vector>& eqn,
161  const word& fieldName
162 ) const
163 {
164  if (debug)
165  {
166  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
167  }
168 
169  if (fieldName == "U")
170  {
171  const surfaceScalarField& rhoPhi =
172  mesh().lookupObject<surfaceScalarField>("rhoPhi");
173 
174  eqn += fvm::Sp(fvc::ddt(rho) + fvc::div(rhoPhi), eqn.psi());
175  }
176 }
177 
178 
180 {}
181 
182 
184 {}
185 
186 
188 {}
189 
190 
192 {
193  return true;
194 }
195 
196 
198 {
199  cavitation_->correct();
200 }
201 
202 
203 // ************************************************************************* //
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:160
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:101
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:79
virtual void correct()
Correct the cavitation model.
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Add implicit/explicit contributions to VoF phase-fraction equation.
Definition: VoFCavitation.C:86
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.
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
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the matrix for implicit and explicit sources.
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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