Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
30#include "fvc.H"
31#include "fvm.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35template<class BasicTurbulenceModel>
37(
38 const word& modelName,
39 const alphaField& alpha,
40 const rhoField& rho,
41 const volVectorField& U,
42 const surfaceScalarField& alphaRhoPhi,
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
64template<class BasicTurbulenceModel>
67 return BasicTurbulenceModel::read();
68}
69
70
71template<class BasicTurbulenceModel>
75 return devRhoReff(this->U_);
76}
77
78
79template<class BasicTurbulenceModel>
82(
83 const volVectorField& U
84) const
85{
87 (
89 (
91 (
92 IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
93 this->runTime_.timeName(),
94 this->mesh_,
97 ),
98 (-(this->alpha_*this->rho_*this->nuEff()))
101 );
102}
103
104
105template<class BasicTurbulenceModel>
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
120template<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
136template<class BasicTurbulenceModel>
138{
139 BasicTurbulenceModel::correct();
140}
141
142
143// ************************************************************************* //
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
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.
virtual bool read()=0
Re-read model coefficients if they have changed.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
A class for managing temporary objects.
Definition tmp.H:75
const surfaceScalarField & alphaRhoPhi_
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
const volScalarField & T
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
volScalarField & alpha