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-2020 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 BasicMomentumTransportModel>
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 )
43 :
45  (
46  type,
47  alpha,
48  rho,
49  U,
50  alphaRhoPhi,
51  phi,
52  transport
53  ),
54 
55  nut_
56  (
57  IOobject
58  (
59  IOobject::groupName("nut", alphaRhoPhi.group()),
60  this->runTime_.timeName(),
61  this->mesh_,
62  IOobject::MUST_READ,
63  IOobject::AUTO_WRITE
64  ),
65  this->mesh_
66  )
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
72 template<class BasicMomentumTransportModel>
74 {
76 }
77 
78 
79 template<class BasicMomentumTransportModel>
82 {
83  tmp<volScalarField> tk(k());
84 
85  // Get list of patchField type names from k
86  wordList patchFieldTypes(tk().boundaryField().types());
87 
88  // For k patchField types which do not have an equivalent for symmTensor
89  // set to calculated
90  forAll(patchFieldTypes, i)
91  {
92  if
93  (
95  ->found(patchFieldTypes[i])
96  )
97  {
99  }
100  }
101 
103  (
104  IOobject::groupName("R", this->alphaRhoPhi_.group()),
105  ((2.0/3.0)*I)*tk() - (nut_)*dev(twoSymm(fvc::grad(this->U_))),
106  patchFieldTypes
107  );
108 }
109 
110 
111 template<class BasicMomentumTransportModel>
113 {
114  correctNut();
115 }
116 
117 
118 template<class BasicMomentumTransportModel>
120 {
122 }
123 
124 
125 // ************************************************************************* //
virtual tmp< volSymmTensorField > sigma() const
Return the Reynolds stress tensor [m^2/s^2].
Definition: eddyViscosity.C:81
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
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
virtual bool read()=0
Re-read model coefficients if they have changed.
Definition: eddyViscosity.C:73
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
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/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
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
phi
Definition: correctPhi.H:3
BasicMomentumTransportModel::rhoField rhoField
Definition: LESModel.H:99
BasicMomentumTransportModel::transportModel transportModel
Definition: LESModel.H:100
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
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
Definition: eddyViscosity.C:34
BasicMomentumTransportModel::alphaField alphaField
Definition: LESModel.H:98
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