cavitationModel.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) 2021 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 "cavitationModel.H"
27 #include "fvcDdt.H"
28 #include "fvcDiv.H"
29 #include "fvmSup.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace twoPhaseChangeModels
36 {
38 }
39 }
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const word& type,
47 )
48 :
49  twoPhaseChangeModel(type, mixture),
50  pSat_("pSat", dimPressure, lookup("pSat"))
51 {}
52 
53 
54 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
55 
58 (
59  volScalarField& alpha
60 ) const
61 {
62  volScalarField::Internal alphalCoeff
63  (
64  1.0/mixture_.rho1()
65  - mixture_.alpha1()()*(1.0/mixture_.rho1() - 1.0/mixture_.rho2())
66  );
67  Pair<tmp<volScalarField>> mDotAlphal = this->mDotAlphal();
68 
69  const volScalarField::Internal vDotcAlphal(alphalCoeff*mDotAlphal[0]());
70  const volScalarField::Internal vDotvAlphal(alphalCoeff*mDotAlphal[1]());
71 
73  (
74  1.0*vDotcAlphal,
75  vDotvAlphal - vDotcAlphal
76  );
77 }
78 
79 
82 (
83  const volScalarField& rho,
84  const volScalarField& gh,
85  volScalarField& p_rgh
86 ) const
87 {
88  dimensionedScalar pCoeff(1.0/mixture_.rho1() - 1.0/mixture_.rho2());
89  Pair<tmp<volScalarField>> mDotP = this->mDotP();
90 
91  const volScalarField vDotcP(pCoeff*mDotP[0]);
92  const volScalarField vDotvP(pCoeff*mDotP[1]);
93 
94  return
95  (vDotvP - vDotcP)*(pSat() - rho*gh)
96  - fvm::Sp(vDotvP - vDotcP, p_rgh);
97 }
98 
99 
102 (
103  const volScalarField& rho,
104  const surfaceScalarField& rhoPhi,
105  volVectorField& U
106 ) const
107 {
108  return fvm::Sp(fvc::ddt(rho) + fvc::div(rhoPhi), U);
109 }
110 
111 
113 {
115  {
116  lookup("pSat") >> pSat_;
117 
118  return true;
119  }
120  else
121  {
122  return false;
123  }
124 }
125 
126 
127 // ************************************************************************* //
defineTypeNameAndDebug(cavitationModel, 0)
Abstract base class for cavitation models.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const dimensionSet dimPressure
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
Calculate the first temporal derivative.
stressControl lookup("compactNormalStress") >> compactNormalStress
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
A class for handling words, derived from string.
Definition: word.H:59
virtual bool read()=0
Read the transportProperties dictionary and update.
virtual Pair< tmp< volScalarField::Internal > > Salpha(volScalarField &alpha) const
Return the cavitation explicit and implicit sources.
Calculate the divergence of the given field.
cavitationModel(const word &type, const immiscibleIncompressibleTwoPhaseMixture &mixture)
Construct for mixture.
virtual tmp< fvVectorMatrix > SU(const volScalarField &rho, const surfaceScalarField &rhoPhi, volVectorField &U) const
Return the cavitation source matrix for the momentum equation.
An immiscible incompressible two-phase mixture transport model.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< fvScalarMatrix > Sp_rgh(const volScalarField &rho, const volScalarField &gh, volScalarField &p_rgh) const
Return the cavitation source matrix.
const volScalarField & gh
A class for managing temporary objects.
Definition: PtrList.H:53
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
virtual bool read()
Read the transportProperties dictionary and update.