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::addAlphaRhoSup
54 (
55  const volScalarField& alpha,
56  const volScalarField& rho,
57  const volScalarField& field,
58  fvMatrix<scalar>& eqn,
59  tmp<volScalarField>
61 ) const
62 {
63  if (debug)
64  {
65  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
66  }
67 
68  const fvMesh& mesh = this->mesh();
69 
70  const phaseSystem::phaseModelPartialList& movingPhases =
71  phase_.fluid().movingPhases();
72 
73  volScalarField::Internal transferRate
74  (
76  (
77  "transferRate",
78  mesh,
80  )
81  );
82  volScalarField::Internal psiTransferRate
83  (
85  (
86  "psiTransferRate",
87  mesh,
88  dimensionedScalar((turbulence_.*psi)()().dimensions()/dimTime, 0)
89  )
90  );
91 
92  forAll(movingPhases, phasei)
93  {
94  if (movingPhases[phasei] != phase_)
95  {
98  (
99  phaseName_
100  );
101 
102  if (notNull(turbulence))
103  {
104  const volScalarField::Internal phaseTransferRate
105  (
106  movingPhases[phasei]
107  *min
108  (
109  turbulence.epsilon()/turbulence.k(),
110  1.0/phase_.time().deltaT()
111  )
112  );
113 
114  transferRate += phaseTransferRate;
115  psiTransferRate += phaseTransferRate*(turbulence.*psi)()();
116  }
117  }
118  }
119 
120  const volScalarField::Internal transferCoeff
121  (
122  max(alphaInversion_ - alpha(), scalar(0))*rho()
123  );
124 
125  eqn += transferCoeff*psiTransferRate;
126  eqn -= fvm::Sp(transferCoeff*transferRate, eqn.psi());
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
131 
133 (
134  const word& sourceName,
135  const word& modelType,
136  const fvMesh& mesh,
137  const dictionary& dict
138 )
139 :
140  fvModel(sourceName, modelType, mesh, dict),
141  phaseName_(dict.lookup("phase")),
142  alphaInversion_("alphaInversion", dimless, dict),
143  phase_
144  (
145  mesh.lookupObject<phaseModel>(IOobject::groupName("alpha", phaseName_))
146  ),
147  turbulence_
148  (
149  mesh.lookupType<phaseCompressible::momentumTransportModel>
150  (
151  phaseName_
152  )
153  )
154 {
155  const word kName(IOobject::groupName("k", phaseName_));
156  const word epsilonName(IOobject::groupName("epsilon", phaseName_));
157  const word omegaName(IOobject::groupName("omega", phaseName_));
158 
159  if (mesh.foundObject<volScalarField>(kName))
160  {
161  fieldNames_.append(kName);
162  }
163 
164  if (mesh.foundObject<volScalarField>(epsilonName))
165  {
166  fieldNames_.append(epsilonName);
167  }
168 
169  if (mesh.foundObject<volScalarField>(omegaName))
170  {
171  fieldNames_.append(omegaName);
172  }
173 }
174 
175 
176 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
177 
179 {
180  return fieldNames_;
181 }
182 
183 
185 (
186  const volScalarField& alpha,
187  const volScalarField& rho,
188  const volScalarField& field,
189  fvMatrix<scalar>& eqn
190 ) const
191 {
192  if (field.name() == IOobject::groupName("k", phaseName_))
193  {
194  addAlphaRhoSup
195  (
196  alpha,
197  rho,
198  field,
199  eqn,
201  );
202  }
203  else if (field.name() == IOobject::groupName("epsilon", phaseName_))
204  {
205  addAlphaRhoSup
206  (
207  alpha,
208  rho,
209  field,
210  eqn,
212  );
213  }
214  else if (field.name() == IOobject::groupName("omega", phaseName_))
215  {
216  addAlphaRhoSup
217  (
218  alpha,
219  rho,
220  field,
221  eqn,
223  );
224  }
225  else
226  {
228  << "Support for field " << field.name() << " is not implemented"
229  << exit(FatalError);
230  }
231 }
232 
233 
235 (
236  const polyTopoChangeMap&
237 )
238 {}
239 
240 
242 {}
243 
244 
246 (
247  const polyDistributionMap&
248 )
249 {}
250 
251 
253 {
254  return true;
255 }
256 
257 
258 // ************************************************************************* //
#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
const word & name() const
Return name.
Definition: IOobject.H:310
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:162
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:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
virtual void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &field, fvMatrix< scalar > &eqn) const
Add contribution to phase k, epsilon or omega equation.
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:163
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:116
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:334
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
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:64
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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))