linearViscousStress.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "linearViscousStress.H"
30 #include "fvc.H"
31 #include "fvm.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class BasicTurbulenceModel>
37 (
38  const word& modelName,
39  const alphaField& alpha,
40  const rhoField& rho,
41  const volVectorField& U,
42  const surfaceScalarField& alphaRhoPhi,
43  const surfaceScalarField& phi,
44  const transportModel& transport,
45  const word& propertiesName
46 )
47 :
48  BasicTurbulenceModel
49  (
50  modelName,
51  alpha,
52  rho,
53  U,
54  alphaRhoPhi,
55  phi,
56  transport,
57  propertiesName
58  )
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
64 template<class BasicTurbulenceModel>
66 {
68 }
69 
70 
71 template<class BasicTurbulenceModel>
74 {
75  return devRhoReff(this->U_);
76 }
77 
78 
79 template<class BasicTurbulenceModel>
82 (
83  const volVectorField& U
84 ) const
85 {
87  (
89  (
90  IOobject
91  (
92  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
93  this->runTime_.timeName(),
94  this->mesh_,
95  IOobject::NO_READ,
96  IOobject::NO_WRITE
97  ),
98  (-(this->alpha_*this->rho_*this->nuEff()))
100  )
101  );
102 }
103 
104 
105 template<class BasicTurbulenceModel>
108 (
110 ) const
111 {
112  return
113  (
114  - fvc::div((this->alpha_*this->rho_*this->nuEff())*dev2(T(fvc::grad(U))))
115  - fvm::laplacian(this->alpha_*this->rho_*this->nuEff(), U)
116  );
117 }
118 
119 
120 template<class BasicTurbulenceModel>
123 (
124  const volScalarField& rho,
126 ) const
127 {
128  return
129  (
130  - fvc::div((this->alpha_*rho*this->nuEff())*dev2(T(fvc::grad(U))))
131  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
132  );
133 }
134 
135 
136 template<class BasicTurbulenceModel>
138 {
140 }
141 
142 
143 // ************************************************************************* //
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:143
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:144
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
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/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
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)
U
Definition: pEqn.H:72
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: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:142
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491