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  (
74  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
75  (-(this->alpha_*this->rho_*this->nuEff()))
76  *dev(twoSymm(fvc::grad(this->U_)))
77  );
78 }
79 
80 
81 template<class BasicTurbulenceModel>
84 (
86 ) const
87 {
88  return
89  (
90  - fvc::div((this->alpha_*this->rho_*this->nuEff())*dev2(T(fvc::grad(U))))
91  - fvm::laplacian(this->alpha_*this->rho_*this->nuEff(), U)
92  );
93 }
94 
95 
96 template<class BasicTurbulenceModel>
99 (
100  const volScalarField& rho,
101  volVectorField& U
102 ) const
103 {
104  return
105  (
106  - fvc::div((this->alpha_*rho*this->nuEff())*dev2(T(fvc::grad(U))))
107  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
108  );
109 }
110 
111 
112 template<class BasicTurbulenceModel>
114 {
116 }
117 
118 
119 // ************************************************************************* //
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
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:89
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:90
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
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(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
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
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:88