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 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  {
39  defineTypeNameAndDebug(phaseTurbulenceStabilisation, 0);
40 
42  (
43  fvModel,
44  phaseTurbulenceStabilisation,
45  dictionary
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  (
99  (
100  momentumTransportModel::typeName,
101  phaseName_
102  )
103  );
104 
105  if (notNull(turbulence))
106  {
107  const volScalarField::Internal phaseTransferRate
108  (
109  movingPhases[phasei]
110  *min
111  (
112  turbulence.epsilon()/turbulence.k(),
113  1.0/phase_.time().deltaT()
114  )
115  );
116 
117  transferRate += phaseTransferRate;
118  psiTransferRate += phaseTransferRate*(turbulence.*psi)()();
119  }
120  }
121  }
122 
123  const volScalarField::Internal transferCoeff
124  (
125  max(alphaInversion_ - alpha(), scalar(0))*rho()
126  );
127 
128  eqn += transferCoeff*psiTransferRate;
129  eqn -= fvm::Sp(transferCoeff*transferRate, eqn.psi());
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
134 
136 (
137  const word& sourceName,
138  const word& modelType,
139  const dictionary& dict,
140  const fvMesh& mesh
141 )
142 :
143  fvModel(sourceName, modelType, dict, mesh),
144  phaseName_(dict.lookup("phase")),
145  alphaInversion_("alphaInversion", dimless, dict),
146  phase_
147  (
148  mesh.lookupObject<phaseModel>(IOobject::groupName("alpha", phaseName_))
149  ),
150  turbulence_
151  (
152  mesh.lookupObject<phaseCompressible::momentumTransportModel>
153  (
154  IOobject::groupName
155  (
156  momentumTransportModel::typeName,
157  phaseName_
158  )
159  )
160  )
161 {
162  const word kName(IOobject::groupName("k", phaseName_));
163  const word epsilonName(IOobject::groupName("epsilon", phaseName_));
164  const word omegaName(IOobject::groupName("omega", phaseName_));
165 
166  if (mesh.foundObject<volScalarField>(kName))
167  {
168  fieldNames_.append(kName);
169  }
170 
171  if (mesh.foundObject<volScalarField>(epsilonName))
172  {
173  fieldNames_.append(epsilonName);
174  }
175 
176  if (mesh.foundObject<volScalarField>(omegaName))
177  {
178  fieldNames_.append(omegaName);
179  }
180 }
181 
182 
183 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
184 
186 {
187  return fieldNames_;
188 }
189 
190 
191 void Foam::fv::phaseTurbulenceStabilisation::addSup
192 (
193  const volScalarField& alpha,
194  const volScalarField& rho,
195  fvMatrix<scalar>& eqn,
196  const word& fieldName
197 ) const
198 {
199  if (fieldName == IOobject::groupName("k", phaseName_))
200  {
201  addSup
202  (
203  alpha,
204  rho,
205  eqn,
207  );
208  }
209  else if (fieldName == IOobject::groupName("epsilon", phaseName_))
210  {
211  addSup
212  (
213  alpha,
214  rho,
215  eqn,
217  );
218  }
219  else if (fieldName == IOobject::groupName("omega", phaseName_))
220  {
221  addSup
222  (
223  alpha,
224  rho,
225  eqn,
227  );
228  }
229  else
230  {
232  << "Support for field " << fieldName << " is not implemented"
233  << exit(FatalError);
234  }
235 }
236 
237 
239 (
240  const polyTopoChangeMap&
241 )
242 {}
243 
244 
245 void Foam::fv::phaseTurbulenceStabilisation::mapMesh(const polyMeshMap& map)
246 {}
247 
248 
250 (
251  const polyDistributionMap&
252 )
253 {}
254 
255 
257 {
258  return true;
259 }
260 
261 
262 // ************************************************************************* //
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
const dimensionSet dimless
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvMesh & mesh
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
compressibleMomentumTransportModel momentumTransportModel
const dimensionSet dimTime
virtual bool movePoints()
Update for mesh motion.
stressControl lookup("compactNormalStress") >> compactNormalStress
Calculate the gradient of the given field.
labelList fv(nPoints)
static word groupName(Name name, const word &group)
const volScalarField & psi
UPtrList< phaseModel > phaseModelPartialList
Definition: phaseSystem.H:84
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:52
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to mixture epsilon or omega equation.
phaseCompressibleMomentumTransportModel momentumTransportModel
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > omega() const =0
Return the turbulence specific dissipation rate.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
phaseTurbulenceStabilisation(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
messageStream Info
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Calculate the matrix for implicit and explicit sources.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Namespace for OpenFOAM.