phaseTurbulenceStabilisation.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 
27 #include "phaseSystem.H"
28 #include "surfaceInterpolate.H"
29 #include "fvcGrad.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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::fv::phaseTurbulenceStabilisation::addSup
54 (
55  const volScalarField& alpha,
56  const volScalarField& rho,
57  fvMatrix<scalar>& eqn,
58  tmp<volScalarField>
60 ) const
61 {
62  if (debug)
63  {
64  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
65  }
66 
67  const fvMesh& mesh = this->mesh();
68 
69  const phaseSystem::phaseModelPartialList& movingPhases =
70  phase_.fluid().movingPhases();
71 
72  volScalarField::Internal transferRate
73  (
75  (
76  "transferRate",
77  mesh,
79  )
80  );
81  volScalarField::Internal psiTransferRate
82  (
84  (
85  "psiTransferRate",
86  mesh,
87  dimensionedScalar((turbulence_.*psi)()().dimensions()/dimTime, 0)
88  )
89  );
90 
91  forAll(movingPhases, phasei)
92  {
93  if (movingPhases[phasei] != phase_)
94  {
97  (
98  phaseName_
99  );
100 
101  if (notNull(turbulence))
102  {
103  const volScalarField::Internal phaseTransferRate
104  (
105  movingPhases[phasei]
106  *min
107  (
108  turbulence.epsilon()/turbulence.k(),
109  1.0/phase_.time().deltaT()
110  )
111  );
112 
113  transferRate += phaseTransferRate;
114  psiTransferRate += phaseTransferRate*(turbulence.*psi)()();
115  }
116  }
117  }
118 
119  const volScalarField::Internal transferCoeff
120  (
121  max(alphaInversion_ - alpha(), scalar(0))*rho()
122  );
123 
124  eqn += transferCoeff*psiTransferRate;
125  eqn -= fvm::Sp(transferCoeff*transferRate, eqn.psi());
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
130 
132 (
133  const word& sourceName,
134  const word& modelType,
135  const fvMesh& mesh,
136  const dictionary& dict
137 )
138 :
139  fvModel(sourceName, modelType, mesh, dict),
140  phaseName_(dict.lookup("phase")),
141  alphaInversion_("alphaInversion", dimless, dict),
142  phase_
143  (
144  mesh.lookupObject<phaseModel>(IOobject::groupName("alpha", phaseName_))
145  ),
146  turbulence_
147  (
148  mesh.lookupType<phaseCompressible::momentumTransportModel>
149  (
150  phaseName_
151  )
152  )
153 {
154  const word kName(IOobject::groupName("k", phaseName_));
155  const word epsilonName(IOobject::groupName("epsilon", phaseName_));
156  const word omegaName(IOobject::groupName("omega", phaseName_));
157 
158  if (mesh.foundObject<volScalarField>(kName))
159  {
160  fieldNames_.append(kName);
161  }
162 
163  if (mesh.foundObject<volScalarField>(epsilonName))
164  {
165  fieldNames_.append(epsilonName);
166  }
167 
168  if (mesh.foundObject<volScalarField>(omegaName))
169  {
170  fieldNames_.append(omegaName);
171  }
172 }
173 
174 
175 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 
178 {
179  return fieldNames_;
180 }
181 
182 
183 void Foam::fv::phaseTurbulenceStabilisation::addSup
184 (
185  const volScalarField& alpha,
186  const volScalarField& rho,
187  fvMatrix<scalar>& eqn,
188  const word& fieldName
189 ) const
190 {
191  if (fieldName == IOobject::groupName("k", phaseName_))
192  {
193  addSup
194  (
195  alpha,
196  rho,
197  eqn,
199  );
200  }
201  else if (fieldName == IOobject::groupName("epsilon", phaseName_))
202  {
203  addSup
204  (
205  alpha,
206  rho,
207  eqn,
209  );
210  }
211  else if (fieldName == IOobject::groupName("omega", phaseName_))
212  {
213  addSup
214  (
215  alpha,
216  rho,
217  eqn,
219  );
220  }
221  else
222  {
224  << "Support for field " << fieldName << " is not implemented"
225  << exit(FatalError);
226  }
227 }
228 
229 
231 (
232  const polyTopoChangeMap&
233 )
234 {}
235 
236 
238 {}
239 
240 
242 (
243  const polyDistributionMap&
244 )
245 {}
246 
247 
249 {
250  return true;
251 }
252 
253 
254 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const Time & time() const
Return time.
Definition: IOobject.C:318
static word groupName(Name name, const word &group)
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:46
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
phaseTurbulenceStabilisation(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volScalarField > omega() const =0
Return the turbulence specific dissipation rate.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
const Type & lookupType(const word &group=word::null) const
Lookup and return the object of the given Type.
bool foundObject(const word &name) const
Is the named Type in registry.
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:127
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:55
UPtrList< phaseModel > phaseModelPartialList
Definition: phaseSystem.H:86
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:306
Calculate the gradient of the given field.
Calculate the matrix for implicit and explicit sources.
const volScalarField & psi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
compressibleMomentumTransportModel momentumTransportModel
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
phaseCompressibleMomentumTransportModel momentumTransportModel
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:58
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))