generalisedNewtonian.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-2021 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 "generalisedNewtonian.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  (
77  (
78  IOobject::modelName("nu", typeName),
79  alphaRhoPhi.group()
80  ),
81  this->runTime_.timeName(),
82  this->mesh_,
85  ),
86  viscosityModel_->nu(this->nu(), strainRate())
87  )
88 {}
89 
90 
91 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
92 
93 template<class BasicMomentumTransportModel>
96 {
97  return sqrt(2.0)*mag(symm(fvc::grad(this->U())));
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
103 template<class BasicMomentumTransportModel>
105 {
106  viscosityModel_->read(this->coeffDict_);
107 
108  return true;
109 }
110 
111 
112 template<class BasicMomentumTransportModel>
115 {
116  return volScalarField::New
117  (
118  IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
119  nu_
120  );
121 }
122 
123 
124 template<class BasicMomentumTransportModel>
127 (
128  const label patchi
129 ) const
130 {
131  return nu_.boundaryField()[patchi];
132 }
133 
134 
135 template<class BasicMomentumTransportModel>
137 {
138  nu_ = viscosityModel_->nu(this->nu(), strainRate());
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 } // End namespace laminarModels
146 } // End namespace Foam
147 
148 // ************************************************************************* //
Foam::surfaceFields.
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 bool read()
Read momentumTransport dictionary.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Calculate the matrix for the laplacian of the field.
Linear viscous stress turbulence model base class.
BasicMomentumTransportModel::transportModel transportModel
generalisedNewtonian(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
BasicMomentumTransportModel::rhoField rhoField
virtual void correct()
Correct the generalisedNewtonian viscosity.
static autoPtr< generalisedNewtonianViscosityModel > New(const dictionary &viscosityProperties)
Return a reference to the selected viscosity model.
Calculate the gradient of the given field.
static word groupName(Name name, const word &group)
BasicMomentumTransportModel::alphaField alphaField
phi
Definition: correctPhi.H:3
Calculate the divergence of the given field.
label patchi
U
Definition: pEqn.H:72
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:247
virtual tmp< volScalarField > strainRate() const
dimensioned< scalar > mag(const dimensioned< Type > &)
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
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.