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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
35 namespace Foam
36 {
37 namespace laminarModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
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  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  (
73  IOobject
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 
93 template<class BasicTurbulenceModel>
95 {
97  {
98  nuM_.readIfPresent(this->coeffDict());
99  lambda_.readIfPresent(this->coeffDict());
100 
101  return true;
102  }
103 
104  return false;
105 }
106 
107 template<class BasicTurbulenceModel>
110 {
111  return sigma_;
112 }
113 
114 template<class BasicTurbulenceModel>
117 {
119  (
121  (
122  IOobject
123  (
124  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
125  this->runTime_.timeName(),
126  this->mesh_,
129  ),
130  this->alpha_*this->rho_*sigma_
131  - (this->alpha_*this->rho_*this->nu())
132  *dev(twoSymm(fvc::grad(this->U_)))
133  )
134  );
135 }
136 
137 
138 template<class BasicTurbulenceModel>
141 (
143 ) const
144 {
145  return
146  (
147  fvc::div
148  (
149  this->alpha_*this->rho_*this->nuM_*fvc::grad(U)
150  )
151  + fvc::div(this->alpha_*this->rho_*sigma_)
152  - fvc::div(this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(U))))
153  - fvm::laplacian(this->alpha_*this->rho_*nu0(), U)
154  );
155 }
156 
157 
158 template<class BasicTurbulenceModel>
161 (
162  const volScalarField& rho,
164 ) const
165 {
166  return
167  (
168  fvc::div
169  (
170  this->alpha_*rho*this->nuM_*fvc::grad(U)
171  )
172  + fvc::div(this->alpha_*rho*sigma_)
173  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
174  - fvm::laplacian(this->alpha_*rho*nu0(), U)
175  );
176 }
177 
178 
179 template<class BasicTurbulenceModel>
181 {
182  // Local references
183  const alphaField& alpha = this->alpha_;
184  const rhoField& rho = this->rho_;
185  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
186  const volVectorField& U = this->U_;
187  volSymmTensorField& sigma = this->sigma_;
188  fv::options& fvOptions(fv::options::New(this->mesh_));
189 
191 
193  const volTensorField& gradU = tgradU();
194 
196  (
197  IOobject
198  (
199  IOobject::groupName("rLambda", this->alphaRhoPhi_.group()),
200  this->runTime_.constant(),
201  this->mesh_
202  ),
203  1.0/(lambda_)
204  );
205 
206  // Note sigma is positive on lhs of momentum eqn
208  (
209  twoSymm(sigma & gradU)
210  - nuM_*rLambda*twoSymm(gradU)
211  );
212 
213  // Viscoelastic stress equation
214  tmp<fvSymmTensorMatrix> sigmaEqn
215  (
217  + fvm::div(alphaRhoPhi, sigma)
218  + fvm::Sp(alpha*rho*rLambda, sigma)
219  ==
220  alpha*rho*P
221  + fvOptions(alpha, rho, sigma)
222  );
223 
224  sigmaEqn.ref().relax();
225  fvOptions.constrain(sigmaEqn.ref());
226  solve(sigmaEqn);
227  fvOptions.correct(sigma_);
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
233 } // End namespace laminarModels
234 } // End namespace Foam
235 
236 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::laminarModels::Maxwell::Maxwell
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:44
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::laminarModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: laminarModel.H:85
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::laminarModels::Maxwell::devRhoReff
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
Definition: Maxwell.C:116
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
Foam::laminarModels::Maxwell::correct
virtual void correct()
Definition: Maxwell.C:180
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::laminarModels::Maxwell::divDevRhoReff
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: Maxwell.C:141
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::dev2
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:117
rho
rho
Definition: readInitialConditions.H:88
Foam::laminarModel
Templated abstract base class for laminar transport models.
Definition: laminarModel.H:52
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::UniformDimensionedField< scalar >
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::dimViscosity
const dimensionSet dimViscosity
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::laminarModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: laminarModel.H:87
Foam::laminarModels::Maxwell::read
virtual bool read()
Read model coefficients if they have changed.
Definition: Maxwell.C:94
uniformDimensionedFields.H
U
U
Definition: pEqn.H:72
Foam::laminarModel::correct
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:333
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Maxwell.H
Foam::laminarModels::Maxwell::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: Maxwell.C:109
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::laminarModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: laminarModel.H:86
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::IOobject::MUST_READ
Definition: IOobject.H:120
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106