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-2018 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
41 (
42  const alphaField& alpha,
43  const rhoField& rho,
44  const volVectorField& U,
45  const surfaceScalarField& alphaRhoPhi,
46  const surfaceScalarField& phi,
47  const transportModel& transport,
48  const word& propertiesName,
49  const word& type
50 )
51 :
53  (
54  type,
55  alpha,
56  rho,
57  U,
58  alphaRhoPhi,
59  phi,
60  transport,
61  propertiesName
62  ),
63 
64  nuM_
65  (
67  (
68  "nuM",
70  this->coeffDict_.lookup("nuM")
71  )
72  ),
73 
74  lambda_
75  (
77  (
78  "lambda",
79  dimTime,
80  this->coeffDict_.lookup("lambda")
81  )
82  ),
83 
84  sigma_
85  (
86  IOobject
87  (
88  IOobject::groupName("sigma", alphaRhoPhi.group()),
89  this->runTime_.timeName(),
90  this->mesh_,
93  ),
94  this->mesh_
95  )
96 {
97  if (type == typeName)
98  {
99  this->printCoeffs(type);
100  }
101 }
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
106 template<class BasicTurbulenceModel>
108 {
110  {
111  nuM_.readIfPresent(this->coeffDict());
112  lambda_.readIfPresent(this->coeffDict());
113 
114  return true;
115  }
116  else
117  {
118  return false;
119  }
120 }
121 
122 template<class BasicTurbulenceModel>
125 {
126  return sigma_;
127 }
128 
129 template<class BasicTurbulenceModel>
132 {
134  (
136  (
137  IOobject
138  (
139  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
140  this->runTime_.timeName(),
141  this->mesh_,
144  ),
145  this->alpha_*this->rho_*sigma_
146  - (this->alpha_*this->rho_*this->nu())
147  *dev(twoSymm(fvc::grad(this->U_)))
148  )
149  );
150 }
151 
152 
153 template<class BasicTurbulenceModel>
156 (
157  volVectorField& U
158 ) const
159 {
160  return
161  (
162  fvc::div
163  (
164  this->alpha_*this->rho_*this->nuM_*fvc::grad(U)
165  )
166  + fvc::div(this->alpha_*this->rho_*sigma_)
167  - fvc::div(this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(U))))
168  - fvm::laplacian(this->alpha_*this->rho_*nu0(), U)
169  );
170 }
171 
172 
173 template<class BasicTurbulenceModel>
176 (
177  const volScalarField& rho,
178  volVectorField& U
179 ) const
180 {
181  return
182  (
183  fvc::div
184  (
185  this->alpha_*rho*this->nuM_*fvc::grad(U)
186  )
187  + fvc::div(this->alpha_*rho*sigma_)
188  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
189  - fvm::laplacian(this->alpha_*rho*nu0(), U)
190  );
191 }
192 
193 
194 template<class BasicTurbulenceModel>
196 {
197  // Local references
198  const alphaField& alpha = this->alpha_;
199  const rhoField& rho = this->rho_;
200  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
201  const volVectorField& U = this->U_;
202  volSymmTensorField& sigma = this->sigma_;
203  fv::options& fvOptions(fv::options::New(this->mesh_));
204 
206 
207  tmp<volTensorField> tgradU(fvc::grad(U));
208  const volTensorField& gradU = tgradU();
209 
211  (
212  IOobject
213  (
214  IOobject::groupName("rLambda", this->alphaRhoPhi_.group()),
215  this->runTime_.constant(),
216  this->mesh_
217  ),
218  1.0/(lambda_)
219  );
220 
221  // Note sigma is positive on lhs of momentum eqn
223  (
224  twoSymm(sigma & gradU)
225  - nuM_*rLambda*twoSymm(gradU)
226  );
227 
228  // Viscoelastic stress equation
229  tmp<fvSymmTensorMatrix> sigmaEqn
230  (
231  fvm::ddt(alpha, rho, sigma)
232  + fvm::div(alphaRhoPhi, sigma)
233  + fvm::Sp(alpha*rho*rLambda, sigma)
234  ==
235  alpha*rho*P
236  + fvOptions(alpha, rho, sigma)
237  );
238 
239  sigmaEqn.ref().relax();
240  fvOptions.constrain(sigmaEqn.ref());
241  solve(sigmaEqn);
242  fvOptions.correct(sigma_);
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 } // End namespace laminarModels
249 } // End namespace Foam
250 
251 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: Maxwell.C:156
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/m2/K4].
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:131
BasicTurbulenceModel::transportModel transportModel
Definition: laminarModel.H:89
BasicTurbulenceModel::rhoField rhoField
Definition: laminarModel.H:88
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
BasicTurbulenceModel::alphaField alphaField
Definition: laminarModel.H:87
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
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:41
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: Maxwell.C:124
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
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:107
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: Maxwell.C:195
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
U
Definition: pEqn.H:72
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:326
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
volScalarField & nu
Namespace for OpenFOAM.