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 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 {
80  return
81  {
82  mixture_.rho1().name(),
83  mixture_.rho2().name()
84  };
85 }
86 
87 
89 (
90  const volScalarField& alpha,
91  const volScalarField& rho,
92  fvMatrix<scalar>& eqn
93 ) const
94 {
95  if (debug)
96  {
97  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
98  }
99 
100  if (&rho == &mixture_.rho1() || &rho == &mixture_.rho2())
101  {
102  const volScalarField& alpha1 = mixture_.alpha1();
103  const volScalarField& alpha2 = mixture_.alpha2();
104 
105  const scalar s = &alpha == &alpha1 ? +1 : -1;
106 
107  // Volume-fraction linearisation
108  if (&alpha == &eqn.psi())
109  {
110  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
111  (
112  cavitation_->mDot12Alpha()
113  );
114  const volScalarField::Internal& mDot1Alpha2 = mDot12Alpha[0]();
115  const volScalarField::Internal& mDot2Alpha1 = mDot12Alpha[1]();
116 
117  eqn +=
118  (&alpha == &alpha1 ? mDot1Alpha2 : mDot2Alpha1)
119  - fvm::Sp(mDot1Alpha2 + mDot2Alpha1, eqn.psi());
120  }
121 
122  // Pressure linearisation
123  else if (eqn.psi().name() == "p_rgh")
124  {
126  (
127  cavitation_->mDot12P()
128  );
129  const volScalarField::Internal& mDot1P = mDot12P[0];
130  const volScalarField::Internal& mDot2P = mDot12P[1];
131 
133  mesh().lookupObject<volScalarField>("rho");
134  const volScalarField::Internal& gh =
135  mesh().lookupObject<volScalarField>("gh");
136 
137  eqn +=
138  fvm::Sp(s*(mDot1P - mDot2P), eqn.psi())
139  + s*(mDot1P - mDot2P)*rho*gh
140  - s*(mDot1P*cavitation_->pSat1() - mDot2P*cavitation_->pSat2());
141  }
142 
143  // Explicit non-linearised value. Used in density predictors and
144  // continuity error terms.
145  else
146  {
147  const Pair<tmp<volScalarField::Internal>> mDot12Alpha
148  (
149  cavitation_->mDot12Alpha()
150  );
151  const volScalarField::Internal mDot1(mDot12Alpha[0]*alpha2);
152  const volScalarField::Internal mDot2(mDot12Alpha[1]*alpha1);
153 
154  eqn += s*(mDot1 - mDot2);
155  }
156  }
157  else
158  {
160  << "Support for field " << alpha.name() << " is not implemented"
161  << exit(FatalError);
162  }
163 }
164 
165 
167 {}
168 
169 
171 {}
172 
173 
175 (
176  const polyDistributionMap&
177 )
178 {}
179 
180 
182 {
183  return true;
184 }
185 
186 
188 {
189  cavitation_->correct();
190 }
191 
192 
193 // ************************************************************************* //
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 rhoFluidThermo-based phases.
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
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:78
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:89
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
#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(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))
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: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