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 | Copyright (C) 2016 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"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace laminarModels
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class BasicTurbulenceModel>
40 (
41  const alphaField& alpha,
42  const rhoField& rho,
43  const volVectorField& U,
44  const surfaceScalarField& alphaRhoPhi,
45  const surfaceScalarField& phi,
46  const transportModel& transport,
47  const word& propertiesName,
48  const word& type
49 )
50 :
52  (
53  type,
54  alpha,
55  rho,
56  U,
57  alphaRhoPhi,
58  phi,
59  transport,
60  propertiesName
61  ),
62 
63  nuM_
64  (
66  (
67  "nuM",
69  this->coeffDict_.lookup("nuM")
70  )
71  ),
72 
73  lambda_
74  (
76  (
77  "lambda",
78  dimTime,
79  this->coeffDict_.lookup("lambda")
80  )
81  ),
82 
83  sigma_
84  (
85  IOobject
86  (
87  IOobject::groupName("sigma", U.group()),
88  this->runTime_.timeName(),
89  this->mesh_,
92  ),
93  this->mesh_
94  )
95 {
96  if (type == typeName)
97  {
98  this->printCoeffs(type);
99  }
100 }
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
105 template<class BasicTurbulenceModel>
107 {
109  {
110  nuM_.readIfPresent(this->coeffDict());
111  lambda_.readIfPresent(this->coeffDict());
112 
113  return true;
114  }
115  else
116  {
117  return false;
118  }
119 }
120 
121 template<class BasicTurbulenceModel>
124 {
125  return sigma_;
126 }
127 
128 template<class BasicTurbulenceModel>
131 {
133  (
135  (
136  IOobject
137  (
138  IOobject::groupName("devRhoReff", this->U_.group()),
139  this->runTime_.timeName(),
140  this->mesh_,
143  ),
144  this->alpha_*this->rho_*sigma_
145  - (this->alpha_*this->rho_*this->nu())
146  *dev(twoSymm(fvc::grad(this->U_)))
147  )
148  );
149 }
150 
151 
152 template<class BasicTurbulenceModel>
155 (
156  volVectorField& U
157 ) const
158 {
159  return
160  (
161  fvc::div
162  (
163  this->alpha_*this->rho_*this->nuM_*fvc::grad(U)
164  )
165  + fvc::div(this->alpha_*this->rho_*sigma_)
166  - fvc::div(this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(U))))
167  - fvm::laplacian(this->alpha_*this->rho_*nu0(), U)
168  );
169 }
170 
171 
172 template<class BasicTurbulenceModel>
175 (
176  const volScalarField& rho,
177  volVectorField& U
178 ) const
179 {
180  return
181  (
182  fvc::div
183  (
184  this->alpha_*rho*this->nuM_*fvc::grad(U)
185  )
186  + fvc::div(this->alpha_*rho*sigma_)
187  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
188  - fvm::laplacian(this->alpha_*rho*nu0(), U)
189  );
190 }
191 
192 
193 template<class BasicTurbulenceModel>
195 {
196  // Local references
197  const alphaField& alpha = this->alpha_;
198  const rhoField& rho = this->rho_;
199  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
200  const volVectorField& U = this->U_;
201  volSymmTensorField& sigma = this->sigma_;
202  fv::options& fvOptions(fv::options::New(this->mesh_));
203 
205 
206  tmp<volTensorField> tgradU(fvc::grad(U));
207  const volTensorField& gradU = tgradU();
208  dimensionedScalar rLambda = 1.0/(lambda_);
209 
210  // Note sigma is positive on lhs of momentum eqn
212  (
213  twoSymm(sigma & gradU)
214  - nuM_*rLambda*twoSymm(gradU)
215  );
216 
217  // Viscoelastic stress equation
218  tmp<fvSymmTensorMatrix> sigmaEqn
219  (
220  fvm::ddt(alpha, rho, sigma)
221  + fvm::div(alphaRhoPhi, sigma)
222  + fvm::Sp(alpha*rho*rLambda, sigma)
223  ==
224  alpha*rho*P
225  + fvOptions(alpha, rho, sigma)
226  );
227 
228  sigmaEqn.ref().relax();
229  fvOptions.constrain(sigmaEqn.ref());
230  solve(sigmaEqn);
231  fvOptions.correct(sigma_);
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace laminarModels
238 } // End namespace Foam
239 
240 // ************************************************************************* //
surfaceScalarField & phi
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:155
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
U
Definition: pEqn.H:83
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:130
fv::options & fvOptions
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:40
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: Maxwell.C:123
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:106
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:485
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: Maxwell.C:194
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
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:103
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.