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-2024 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 "fvcGrad.H"
28 #include "fvcDiv.H"
29 #include "fvmLaplacian.H"
30 #include "fvcSnGrad.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class BasicMomentumTransportModel>
36 (
37  const word& modelName,
38  const alphaField& alpha,
39  const rhoField& rho,
40  const volVectorField& U,
41  const surfaceScalarField& alphaRhoPhi,
42  const surfaceScalarField& phi,
43  const viscosity& viscosity
44 )
45 :
46  BasicMomentumTransportModel
47  (
48  modelName,
49  alpha,
50  rho,
51  U,
52  alphaRhoPhi,
53  phi,
54  viscosity
55  )
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
61 template<class BasicMomentumTransportModel>
63 {
65 }
66 
67 
68 template<class BasicMomentumTransportModel>
71 {
72  const surfaceScalarField alphaRhoNuEff
73  (
74  fvc::interpolate(this->alpha_*this->rho_*this->nuEff())
75  );
76 
78  (
79  this->groupName("devTau"),
80  -alphaRhoNuEff
81  *(
82  fvc::dotInterpolate(this->mesh().nf(), dev2(T(fvc::grad(this->U_))))
83  + fvc::snGrad(this->U_)
84  )
85  );
86 }
87 
88 
89 template<class BasicMomentumTransportModel>
90 template<class RhoFieldType>
93 (
94  const RhoFieldType& rho,
96 ) const
97 {
98  const surfaceScalarField alphaRhoNuEff
99  (
100  fvc::interpolate(this->alpha_*rho*this->nuEff())
101  );
102 
103  const fvVectorMatrix divDevTauCorr
104  (
105  this->divDevTauCorr
106  (
107  -alphaRhoNuEff
108  *fvc::dotInterpolate(this->mesh().Sf(), dev2(T(fvc::grad(U)))),
110  )
111  );
112 
113  return divDevTauCorr - fvm::laplacian(alphaRhoNuEff, U);
114 }
115 
116 
117 template<class BasicMomentumTransportModel>
120 (
122 ) const
123 {
124  return DivDevTau(this->rho_, U);
125 }
126 
127 
128 template<class BasicMomentumTransportModel>
131 (
132  const volScalarField& rho,
134 ) const
135 {
136  return DivDevTau(rho, U);
137 }
138 
139 
140 template<class BasicMomentumTransportModel>
142 {
144 }
145 
146 
147 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
BasicMomentumTransportModel::alphaField alphaField
tmp< fvVectorMatrix > DivDevTau(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
linearViscousStress(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual tmp< surfaceVectorField > devTau() const
Return the effective surface stress.
virtual bool read()=0
Re-read model coefficients if they have changed.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const scalar nuEff
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the snGrad of the given volField.
Calculate the matrix for the laplacian of the field.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate(const surfaceVectorField &Sf, const VolField< Type > &tvf)
Interpolate field onto faces.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
void dev2(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)