Loading...
Searching...
No Matches
generalizedNewtonian.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) 2018-2020 OpenFOAM Foundation
9 Copyright (C) 2020-2021 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 "volFields.H"
31#include "surfaceFields.H"
32#include "fvcGrad.H"
33#include "fvcDiv.H"
34#include "fvmLaplacian.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
40namespace laminarModels
41{
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45template<class BasicMomentumTransportModel>
47(
48 const alphaField& alpha,
49 const rhoField& rho,
50 const volVectorField& U,
51 const surfaceScalarField& alphaRhoPhi,
53 const transportModel& transport,
54 const word& propertiesName
55)
56:
58 (
60 alpha,
61 rho,
62 U,
63 alphaRhoPhi,
64 phi,
65 transport,
66 propertiesName
67 ),
68
70 (
72 (
73 this->coeffDict_
74 )
75 ),
76
77 nu_
78 (
80 (
81 IOobject::groupName("generalizedNewtonian:nu", alphaRhoPhi.group()),
82 this->runTime_.timeName(),
83 this->mesh_,
86 ),
87 viscosityModel_->nu(this->nu(), strainRate())
88 )
89{}
90
91
92// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
93
94template<class BasicMomentumTransportModel>
97{
98 return sqrt(scalar{2})*mag(symm(fvc::grad(this->U())));
99}
100
101
102// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103
104template<class BasicMomentumTransportModel>
106{
107 viscosityModel_->read(this->coeffDict_);
109 return true;
110}
111
112
113template<class BasicMomentumTransportModel>
116{
118 (
119 IOobject::groupName("nut", this->alphaRhoPhi_.group()),
120 this->mesh_,
122 );
123}
124
125
126template<class BasicMomentumTransportModel>
129(
130 const label patchi
131) const
133 return tmp<scalarField>::New(this->mesh_.boundary()[patchi].size(), Zero);
134}
135
136
137template<class BasicMomentumTransportModel>
140{
142 (
143 IOobject::groupName("nuEff", this->alphaRhoPhi_.group()),
145 );
146}
147
148
149template<class BasicMomentumTransportModel>
152(
153 const label patchi
154) const
155{
156 return nu_.boundaryField()[patchi];
157}
158
159
160template<class BasicMomentumTransportModel>
162{
163 nu_ = viscosityModel_->nu(this->nu(), strainRate());
165}
166
167
168// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169
170} // End namespace laminarModels
171} // End namespace Foam
172
173// ************************************************************************* //
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word group(const word &name)
Return group (extension part of name).
Definition IOobject.C:300
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual void correct()
Correct the laminar transport.
An abstract base class for generalized Newtonian viscosity models.
static autoPtr< generalizedNewtonianViscosityModel > New(const dictionary &dict)
Select a viscosity model.
BasicMomentumTransportModel::alphaField alphaField
generalizedNewtonian(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Construct from components.
autoPtr< generalizedNewtonianViscosityModel > viscosityModel_
Run-time selectable non-Newtonian viscosity model.
virtual void correct()
Correct the generalizedNewtonian viscosity.
volScalarField nu_
The non-Newtonian viscosity field.
virtual tmp< volScalarField > strainRate() const
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
BasicMomentumTransportModel::transportModel transportModel
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity,.
virtual bool read()
Read turbulence (momentumTransport) dictionary.
BasicMomentumTransportModel::rhoField rhoField
static autoPtr< generalizedNewtonian > New(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Return a reference to the selected turbulence model.
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)
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the matrix for the laplacian of the field.
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
Namespace for OpenFOAM.
const dimensionSet dimViscosity
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & nu
volScalarField & alpha
Foam::surfaceFields.