linearViscousStress.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 "linearViscousStress.H"
27 #include "fvc.H"
28 #include "fvm.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class BasicTurbulenceModel>
34 (
35  const word& modelName,
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 :
45  BasicTurbulenceModel
46  (
47  modelName,
48  alpha,
49  rho,
50  U,
51  alphaRhoPhi,
52  phi,
53  transport,
54  propertiesName
55  )
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
61 template<class BasicTurbulenceModel>
63 {
65 }
66 
67 
68 template<class BasicTurbulenceModel>
71 {
73  (
75  (
76  IOobject
77  (
78  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
79  this->runTime_.timeName(),
80  this->mesh_,
81  IOobject::NO_READ,
82  IOobject::NO_WRITE
83  ),
84  (-(this->alpha_*this->rho_*this->nuEff()))
85  *dev(twoSymm(fvc::grad(this->U_)))
86  )
87  );
88 }
89 
90 
91 template<class BasicTurbulenceModel>
94 (
96 ) const
97 {
98  return
99  (
100  - fvc::div((this->alpha_*this->rho_*this->nuEff())*dev2(T(fvc::grad(U))))
101  - fvm::laplacian(this->alpha_*this->rho_*this->nuEff(), U)
102  );
103 }
104 
105 
106 template<class BasicTurbulenceModel>
109 (
110  const volScalarField& rho,
111  volVectorField& U
112 ) const
113 {
114  return
115  (
116  - fvc::div((this->alpha_*rho*this->nuEff())*dev2(T(fvc::grad(U))))
117  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
118  );
119 }
120 
121 
122 template<class BasicTurbulenceModel>
124 {
126 }
127 
128 
129 // ************************************************************************* //
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
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(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/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:72
virtual bool read()=0
Re-read model coefficients if they have changed.
const volScalarField & T
linearViscousStress(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.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
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
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99