Loading...
Searching...
No Matches
Maxwell.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) 2016-2017 OpenFOAM Foundation
9 Copyright (C) 2019-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
29#include "Maxwell.H"
30#include "fvOptions.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
37namespace laminarModels
38{
39
40// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41
42template<class BasicTurbulenceModel>
44(
45 const alphaField& alpha,
46 const rhoField& rho,
47 const volVectorField& U,
48 const surfaceScalarField& alphaRhoPhi,
50 const transportModel& transport,
51 const word& propertiesName,
52 const word& type
53)
54:
56 (
57 type,
58 alpha,
59 rho,
60 U,
61 alphaRhoPhi,
62 phi,
63 transport,
64 propertiesName
65 ),
66
67 nuM_("nuM", dimViscosity, this->coeffDict_),
68
69 lambda_("lambda", dimTime, this->coeffDict_),
70
71 sigma_
72 (
74 (
75 IOobject::groupName("sigma", alphaRhoPhi.group()),
76 this->runTime_.timeName(),
77 this->mesh_,
80 ),
81 this->mesh_
82 )
83{
84 if (type == typeName)
85 {
86 this->printCoeffs(type);
87 }
88}
89
90
91// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92
93template<class BasicTurbulenceModel>
95{
97 {
98 nuM_.readIfPresent(this->coeffDict());
99 lambda_.readIfPresent(this->coeffDict());
100
101 return true;
102 }
104 return false;
105}
106
107
108template<class BasicTurbulenceModel>
112 return sigma_;
113}
114
115
116template<class BasicTurbulenceModel>
120 return devRhoReff(this->U_);
121}
122
123
124template<class BasicTurbulenceModel>
127{
129 (
131 (
133 (
134 IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
135 this->runTime_.timeName(),
136 this->mesh_,
139 ),
140 this->alpha_*this->rho_*sigma_
141 - (this->alpha_*this->rho_*this->nu())
144 );
145}
146
147
148template<class BasicTurbulenceModel>
151(
153) const
154{
155 return
156 (
158 (
159 this->alpha_*this->rho_*this->nuM_*fvc::grad(U)
160 )
161 + fvc::div(this->alpha_*this->rho_*sigma_)
162 - fvc::div(this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(U))))
163 - fvm::laplacian(this->alpha_*this->rho_*nu0(), U)
164 );
165}
166
167
168template<class BasicTurbulenceModel>
171(
172 const volScalarField& rho,
174) const
175{
176 return
177 (
179 (
180 this->alpha_*rho*this->nuM_*fvc::grad(U)
181 )
182 + fvc::div(this->alpha_*rho*sigma_)
183 - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
184 - fvm::laplacian(this->alpha_*rho*nu0(), U)
185 );
186}
187
188
189template<class BasicTurbulenceModel>
191{
192 // Local references
193 const alphaField& alpha = this->alpha_;
194 const rhoField& rho = this->rho_;
195 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
196 const volVectorField& U = this->U_;
197 volSymmTensorField& sigma = this->sigma_;
199
201
203 const volTensorField& gradU = tgradU();
204
206 (
208 (
209 IOobject::groupName("rLambda", this->alphaRhoPhi_.group()),
210 this->runTime_.constant(),
211 this->mesh_
212 ),
213 1.0/(lambda_)
214 );
215
216 // Note sigma is positive on lhs of momentum eqn
218 (
219 twoSymm(sigma & gradU)
220 - nuM_*rLambda*twoSymm(gradU)
221 );
222
223 // Viscoelastic stress equation
225 (
227 + fvm::div(alphaRhoPhi, sigma)
228 + fvm::Sp(alpha*rho*rLambda, sigma)
229 ==
230 alpha*rho*P
232 );
233
234 sigmaEqn.ref().relax();
235 fvOptions.constrain(sigmaEqn.ref());
236 solve(sigmaEqn);
237 fvOptions.correct(sigma_);
238}
239
240
241// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242
243} // End namespace laminarModels
244} // End namespace Foam
245
246// ************************************************************************* //
fv::options & fvOptions
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ 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.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
dictionary coeffDict_
Model coefficients dictionary.
virtual void printCoeffs(const word &type)
Print model coefficients.
virtual void correct()
Correct the laminar transport.
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
laminarModel(const laminarModel &)=delete
No copy construct.
virtual bool read()
Read model coefficients if they have changed.
volSymmTensorField sigma_
Definition Maxwell.H:87
BasicTurbulenceModel::alphaField alphaField
Definition Maxwell.H:103
BasicTurbulenceModel::rhoField rhoField
Definition Maxwell.H:104
dimensionedScalar lambda_
Definition Maxwell.H:82
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and related properties.
Definition Maxwell.C:183
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition Maxwell.C:103
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
Definition Maxwell.C:111
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition Maxwell.C:144
BasicTurbulenceModel::transportModel transportModel
Definition Maxwell.H:105
dimensionedScalar nuM_
Definition Maxwell.H:81
tmp< volScalarField > nu0() const
Return the turbulence viscosity.
Definition Maxwell.H:95
Maxwell(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition Maxwell.C:37
virtual bool read()
Read model coefficients if they have changed.
Definition Maxwell.C:87
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
const volScalarField & T
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
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)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
const dimensionSet dimViscosity
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UniformDimensionedField< scalar > uniformDimensionedScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
volScalarField & nu
volScalarField & alpha
CEqn solve()
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Various UniformDimensionedField types.