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 "fvmSup.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace fv
36  {
37  namespace compressible
38  {
40 
42  (
43  fvModel,
46  );
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<compressibleTwoPhaseVoFMixture>
67  (
68  "phaseProperties"
69  )
70  ),
71 
72  cavitation_(Foam::compressible::cavitationModel::New(dict, mixture_)),
73 
74  alphaName_(mixture_.alpha1().name())
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  return {alphaName_, "p_rgh"};
83 }
84 
85 
87 (
88  fvMatrix<scalar>& eqn,
89  const word& fieldName
90 ) const
91 {
92  if (debug)
93  {
94  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
95  }
96 
97  if (fieldName == alphaName_)
98  {
99  const volScalarField::Internal alpha1Coeff
100  (
101  1.0/mixture_.rho1()()
102  - mixture_.alpha1()()*(1.0/mixture_.rho1()() - 1.0/mixture_.rho2()())
103  );
104 
105  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
106  (
107  cavitation_->mDot12Alpha()
108  );
109 
110  const volScalarField::Internal vDot1Alpha(alpha1Coeff*mDot12Alpha[0]());
111  const volScalarField::Internal vDot2Alpha(alpha1Coeff*mDot12Alpha[1]());
112 
113  eqn += fvm::Sp(-vDot2Alpha - vDot1Alpha, eqn.psi()) + vDot1Alpha;
114  }
115 }
116 
117 
119 (
120  const volScalarField&,
121  fvMatrix<scalar>& eqn,
122  const word& fieldName
123 ) const
124 {
125  if (debug)
126  {
127  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
128  }
129 
130  if (fieldName == "p_rgh")
131  {
133  mesh().lookupObject<volScalarField>("rho");
134 
135  const volScalarField::Internal& gh =
136  mesh().lookupObject<volScalarField>("gh");
137 
138  const volScalarField::Internal pCoeff
139  (
140  1.0/mixture_.rho1()() - 1.0/mixture_.rho2()()
141  );
142 
144  (
145  cavitation_->mDot12P()
146  );
147 
148  const volScalarField::Internal vDot1P(pCoeff*mDot12P[0]);
149  const volScalarField::Internal vDot2P(pCoeff*mDot12P[1]);
150 
151  eqn +=
152  vDot2P*cavitation_->pSat1() - vDot1P*cavitation_->pSat2()
153  - (vDot2P - vDot1P)*rho*gh
154  - fvm::Sp(vDot2P - vDot1P, eqn.psi());
155  }
156 }
157 
158 
160 {}
161 
162 
164 {}
165 
166 
168 (
169  const polyDistributionMap&
170 )
171 {}
172 
173 
175 {
176  return true;
177 }
178 
179 
181 {
182  cavitation_->correct();
183 }
184 
185 
186 // ************************************************************************* //
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.
Class to represent a mixture of two rhoThermo-based phases.
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
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:80
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:87
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:55
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 matrix for implicit and explicit sources.
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.
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