eddyViscosity.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) 2013-2018 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 "eddyViscosity.H"
27 #include "fvc.H"
28 #include "fvm.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class BasicTurbulenceModel>
34 (
35  const word& type,
36  const alphaField& alpha,
37  const rhoField& rho,
38  const volVectorField& U,
39  const surfaceScalarField& alphaRhoPhi,
40  const surfaceScalarField& phi,
41  const transportModel& transport,
42  const word& propertiesName
43 )
44 :
46  (
47  type,
48  alpha,
49  rho,
50  U,
51  alphaRhoPhi,
52  phi,
53  transport,
54  propertiesName
55  ),
56 
57  nut_
58  (
59  IOobject
60  (
61  IOobject::groupName("nut", alphaRhoPhi.group()),
62  this->runTime_.timeName(),
63  this->mesh_,
64  IOobject::MUST_READ,
65  IOobject::AUTO_WRITE
66  ),
67  this->mesh_
68  )
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
74 template<class BasicTurbulenceModel>
76 {
78 }
79 
80 
81 template<class BasicTurbulenceModel>
84 {
85  tmp<volScalarField> tk(k());
86 
87  // Get list of patchField type names from k
88  wordList patchFieldTypes(tk().boundaryField().types());
89 
90  // For k patchField types which do not have an equivalent for symmTensor
91  // set to calculated
92  forAll(patchFieldTypes, i)
93  {
94  if
95  (
97  ->found(patchFieldTypes[i])
98  )
99  {
101  }
102  }
103 
105  (
106  IOobject::groupName("R", this->alphaRhoPhi_.group()),
107  ((2.0/3.0)*I)*tk() - (nut_)*dev(twoSymm(fvc::grad(this->U_))),
108  patchFieldTypes
109  );
110 }
111 
112 
113 template<class BasicTurbulenceModel>
115 {
116  correctNut();
117 }
118 
119 
120 template<class BasicTurbulenceModel>
122 {
124 }
125 
126 
127 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
surfaceScalarField & phi
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:89
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
Linear viscous stress turbulence model base class.
label k
Boltzmann constant.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: eddyViscosity.C:83
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:90
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
virtual bool read()=0
Re-read model coefficients if they have changed.
Definition: eddyViscosity.C:75
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:93
A class for handling words, derived from string.
Definition: word.H:59
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-7/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
U
Definition: pEqn.H:72
virtual void validate()
Validate the turbulence fields after construction.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
bool found
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Construct from components.
Definition: eddyViscosity.C:34
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:88