Maxwell.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) 2016-2019 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 "Maxwell.H"
27 #include "fvOptions.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace laminarModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
41 {
43  (
45  (
46  sigma_,
47  dimVolume*this->rho_.dimensions()*sigma_.dimensions()/dimTime
48  )
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 template<class BasicTurbulenceModel>
57 (
58  const alphaField& alpha,
59  const rhoField& rho,
60  const volVectorField& U,
61  const surfaceScalarField& alphaRhoPhi,
62  const surfaceScalarField& phi,
63  const transportModel& transport,
64  const word& propertiesName,
65  const word& type
66 )
67 :
69  (
70  type,
71  alpha,
72  rho,
73  U,
74  alphaRhoPhi,
75  phi,
76  transport,
77  propertiesName
78  ),
79 
80  nuM_
81  (
83  (
84  "nuM",
86  this->coeffDict_.lookup("nuM")
87  )
88  ),
89 
90  lambda_
91  (
93  (
94  "lambda",
95  dimTime,
96  this->coeffDict_.lookup("lambda")
97  )
98  ),
99 
100  sigma_
101  (
102  IOobject
103  (
104  IOobject::groupName("sigma", alphaRhoPhi.group()),
105  this->runTime_.timeName(),
106  this->mesh_,
109  ),
110  this->mesh_
111  )
112 {
113  if (type == typeName)
114  {
115  this->printCoeffs(type);
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
122 template<class BasicTurbulenceModel>
124 {
126  {
127  nuM_.readIfPresent(this->coeffDict());
128  lambda_.readIfPresent(this->coeffDict());
129 
130  return true;
131  }
132  else
133  {
134  return false;
135  }
136 }
137 
138 template<class BasicTurbulenceModel>
141 {
142  return sigma_;
143 }
144 
145 template<class BasicTurbulenceModel>
148 {
150  (
151  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
152  this->alpha_*this->rho_*sigma_
153  - (this->alpha_*this->rho_*this->nu())
154  *dev(twoSymm(fvc::grad(this->U_)))
155  );
156 }
157 
158 
159 template<class BasicTurbulenceModel>
162 (
163  volVectorField& U
164 ) const
165 {
166  return
167  (
168  fvc::div
169  (
170  this->alpha_*this->rho_*this->nuM_*fvc::grad(U)
171  )
172  + fvc::div(this->alpha_*this->rho_*sigma_)
173  - fvc::div(this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(U))))
174  - fvm::laplacian(this->alpha_*this->rho_*nu0(), U)
175  );
176 }
177 
178 
179 template<class BasicTurbulenceModel>
182 (
183  const volScalarField& rho,
184  volVectorField& U
185 ) const
186 {
187  return
188  (
189  fvc::div
190  (
191  this->alpha_*rho*this->nuM_*fvc::grad(U)
192  )
193  + fvc::div(this->alpha_*rho*sigma_)
194  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
195  - fvm::laplacian(this->alpha_*rho*nu0(), U)
196  );
197 }
198 
199 
200 template<class BasicTurbulenceModel>
202 {
203  // Local references
204  const alphaField& alpha = this->alpha_;
205  const rhoField& rho = this->rho_;
206  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
207  const volVectorField& U = this->U_;
208  volSymmTensorField& sigma = this->sigma_;
209  fv::options& fvOptions(fv::options::New(this->mesh_));
210 
212 
213  tmp<volTensorField> tgradU(fvc::grad(U));
214  const volTensorField& gradU = tgradU();
215 
217  (
218  IOobject
219  (
220  IOobject::groupName("rLambda", this->alphaRhoPhi_.group()),
221  this->runTime_.constant(),
222  this->mesh_
223  ),
224  1.0/(lambda_)
225  );
226 
227  // Note sigma is positive on lhs of momentum eqn
229  (
230  twoSymm(sigma & gradU)
231  - nuM_*rLambda*twoSymm(gradU)
232  );
233 
234  // Viscoelastic stress equation
235  tmp<fvSymmTensorMatrix> sigmaEqn
236  (
237  fvm::ddt(alpha, rho, sigma)
238  + fvm::div(alphaRhoPhi, sigma)
239  + fvm::Sp(alpha*rho*rLambda, sigma)
240  ==
241  alpha*rho*P
242  + sigmaSource()
243  + fvOptions(alpha, rho, sigma)
244  );
245 
246  sigmaEqn.ref().relax();
247  fvOptions.constrain(sigmaEqn.ref());
248  solve(sigmaEqn);
249  fvOptions.correct(sigma_);
250 }
251 
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 } // End namespace laminarModels
256 } // End namespace Foam
257 
258 // ************************************************************************* //
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: Maxwell.C:162
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Templated abstract base class for laminar transport models.
Definition: laminarModel.H:49
fv::options & fvOptions
surfaceScalarField & phi
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const dimensionSet dimViscosity
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
fvMatrix< symmTensor > fvSymmTensorMatrix
Definition: fvMatricesFwd.H:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Finite-volume options.
Definition: fvOptions.H:52
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
Definition: Maxwell.C:147
BasicTurbulenceModel::transportModel transportModel
Definition: laminarModel.H:78
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
BasicTurbulenceModel::rhoField rhoField
Definition: laminarModel.H:77
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
BasicTurbulenceModel::alphaField alphaField
Definition: laminarModel.H:76
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
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:57
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: Maxwell.C:140
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
virtual tmp< fvSymmTensorMatrix > sigmaSource() const
Definition: Maxwell.C:40
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
virtual bool read()
Read model coefficients if they have changed.
Definition: Maxwell.C:123
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: Maxwell.C:201
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
U
Definition: pEqn.H:72
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:274
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
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
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
volScalarField & nu
Namespace for OpenFOAM.