generalizedNewtonian.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) 2018-2020 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 "generalizedNewtonian.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "fvcGrad.H"
30 #include "fvcDiv.H"
31 #include "fvmLaplacian.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace laminarModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicMomentumTransportModel>
44 (
45  const alphaField& alpha,
46  const rhoField& rho,
47  const volVectorField& U,
48  const surfaceScalarField& alphaRhoPhi,
49  const surfaceScalarField& phi,
50  const transportModel& transport
51 )
52 :
54  (
55  typeName,
56  alpha,
57  rho,
58  U,
59  alphaRhoPhi,
60  phi,
61  transport
62  ),
63 
64  viscosityModel_
65  (
67  (
68  this->coeffDict_
69  )
70  ),
71 
72  nu_
73  (
74  IOobject
75  (
76  IOobject::groupName("generalizedNewtonian:nu", alphaRhoPhi.group()),
77  this->runTime_.timeName(),
78  this->mesh_,
81  ),
82  viscosityModel_->nu(this->nu(), strainRate())
83  )
84 {}
85 
86 
87 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
88 
89 template<class BasicMomentumTransportModel>
92 {
93  return sqrt(2.0)*mag(symm(fvc::grad(this->U())));
94 }
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class BasicMomentumTransportModel>
101 {
102  viscosityModel_->read(this->coeffDict_);
103 
104  return true;
105 }
106 
107 
108 template<class BasicMomentumTransportModel>
111 {
112  return volScalarField::New
113  (
114  IOobject::groupName("nut", this->alphaRhoPhi_.group()),
115  this->mesh_,
117  );
118 }
119 
120 
121 template<class BasicMomentumTransportModel>
124 (
125  const label patchi
126 ) const
127 {
128  return tmp<scalarField>
129  (
130  new scalarField(this->mesh_.boundary()[patchi].size(), 0.0)
131  );
132 }
133 
134 
135 template<class BasicMomentumTransportModel>
138 {
139  return volScalarField::New
140  (
141  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
142  nu_
143  );
144 }
145 
146 
147 template<class BasicMomentumTransportModel>
150 (
151  const label patchi
152 ) const
153 {
154  return nu_.boundaryField()[patchi];
155 }
156 
157 
158 template<class BasicMomentumTransportModel>
161 {
162  return volScalarField::New
163  (
164  IOobject::groupName("k", this->alphaRhoPhi_.group()),
165  this->mesh_,
166  dimensionedScalar(sqr(this->U_.dimensions()), 0)
167  );
168 }
169 
170 
171 template<class BasicMomentumTransportModel>
174 {
175  return volScalarField::New
176  (
177  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
178  this->mesh_,
179  dimensionedScalar(sqr(this->U_.dimensions())/dimTime, 0)
180  );
181 }
182 
183 
184 template<class BasicMomentumTransportModel>
187 {
189  (
190  IOobject::groupName("R", this->alphaRhoPhi_.group()),
191  this->mesh_,
192  dimensionedSymmTensor(sqr(this->U_.dimensions()), Zero)
193  );
194 }
195 
196 
197 template<class BasicMomentumTransportModel>
199 {
200  nu_ = viscosityModel_->nu(this->nu(), strainRate());
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace laminarModels
208 } // End namespace Foam
209 
210 // ************************************************************************* //
Foam::surfaceFields.
static autoPtr< generalizedNewtonianViscosityModel > New(const dictionary &viscosityProperties)
Return a reference to the selected viscosity model.
BasicMomentumTransportModel::transportModel transportModel
BasicMomentumTransportModel::alphaField alphaField
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual void correct()
Correct the generalizedNewtonian viscosity.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
const dimensionSet dimViscosity
dimensionedScalar sqrt(const dimensionedScalar &ds)
Calculate the matrix for the laplacian of the field.
Linear viscous stress turbulence model base class.
generalizedNewtonian(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
phi
Definition: pEqn.H:104
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate,.
BasicMomentumTransportModel::rhoField rhoField
Calculate the gradient of the given field.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
static const zero Zero
Definition: zero.H:97
virtual bool read()
Read momentumTransport dictionary.
Calculate the divergence of the given field.
virtual tmp< volScalarField > strainRate() const
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity,.
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
label patchi
U
Definition: pEqn.H:72
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:271
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
Namespace for OpenFOAM.