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-2020 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 BasicMomentumTransportModel>
42 (
43  const word& name,
44  const dimensionSet& dims
45 ) const
46 {
47  PtrList<dimensionedScalar> modeCoeffs(nModes_);
48 
49  if (modeCoefficients_.size())
50  {
51  if (this->coeffDict().found(name))
52  {
53  IOWarningInFunction(this->coeffDict())
54  << "Using 'modes' list, '" << name << "' entry will be ignored."
55  << endl;
56  }
57 
58  forAll(modeCoefficients_, modei)
59  {
60  modeCoeffs.set
61  (
62  modei,
64  (
65  name,
66  dims,
67  modeCoefficients_[modei].lookup(name)
68  )
69  );
70  }
71  }
72  else
73  {
74  modeCoeffs.set
75  (
76  0,
78  (
79  name,
80  dims,
81  this->coeffDict_.lookup(name)
82  )
83  );
84  }
85 
86  return modeCoeffs;
87 }
88 
89 
90 template<class BasicMomentumTransportModel>
92 (
93  const label modei,
94  volSymmTensorField& sigma
95 ) const
96 {
97  return -fvm::Sp(this->alpha_*this->rho_/lambdas_[modei], sigma);
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
103 template<class BasicMomentumTransportModel>
105 (
106  const alphaField& alpha,
107  const rhoField& rho,
108  const volVectorField& U,
109  const surfaceScalarField& alphaRhoPhi,
110  const surfaceScalarField& phi,
111  const transportModel& transport,
112  const word& type
113 )
114 :
116  (
117  type,
118  alpha,
119  rho,
120  U,
121  alphaRhoPhi,
122  phi,
123  transport
124  ),
125 
126  modeCoefficients_
127  (
128  this->coeffDict().found("modes")
130  (
131  this->coeffDict().lookup("modes")
132  )
134  ),
135 
136  nModes_(modeCoefficients_.size() ? modeCoefficients_.size() : 1),
137 
138  nuM_
139  (
141  (
142  "nuM",
143  dimViscosity,
144  this->coeffDict_.lookup("nuM")
145  )
146  ),
147 
148  lambdas_(readModeCoefficients("lambda", dimTime)),
149 
150  sigma_
151  (
152  IOobject
153  (
154  IOobject::groupName("sigma", alphaRhoPhi.group()),
155  this->runTime_.timeName(),
156  this->mesh_,
159  ),
160  this->mesh_
161  )
162 {
163  if (nModes_ > 1)
164  {
165  sigmas_.setSize(nModes_);
166 
167  forAll(sigmas_, modei)
168  {
169  IOobject header
170  (
171  IOobject::groupName("sigma" + name(modei), alphaRhoPhi.group()),
172  this->runTime_.timeName(),
173  this->mesh_,
175  );
176 
177  // Check if mode field exists and can be read
178  if (header.typeHeaderOk<volSymmTensorField>(true))
179  {
180  Info<< " Reading mode stress field "
181  << header.name() << endl;
182 
183  sigmas_.set
184  (
185  modei,
187  (
188  IOobject
189  (
190  header.name(),
191  this->runTime_.timeName(),
192  this->mesh_,
195  ),
196  this->mesh_
197  )
198  );
199  }
200  else
201  {
202  sigmas_.set
203  (
204  modei,
206  (
207  IOobject
208  (
209  header.name(),
210  this->runTime_.timeName(),
211  this->mesh_,
214  ),
215  sigma_
216  )
217  );
218  }
219  }
220  }
221 
222  if (type == typeName)
223  {
224  this->printCoeffs(type);
225  }
226 }
227 
228 
229 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
230 
231 template<class BasicMomentumTransportModel>
233 {
235  {
236  if (modeCoefficients_.size())
237  {
238  this->coeffDict().lookup("modes") >> modeCoefficients_;
239  }
240 
241  nuM_.readIfPresent(this->coeffDict());
242 
243  lambdas_ = readModeCoefficients("lambda", dimTime);
244 
245  return true;
246  }
247  else
248  {
249  return false;
250  }
251 }
252 
253 
254 template<class BasicMomentumTransportModel>
257 {
258  return sigma_;
259 }
260 
261 
262 template<class BasicMomentumTransportModel>
265 {
267  (
268  IOobject::groupName("devTau", this->alphaRhoPhi_.group()),
269  this->alpha_*this->rho_*sigma_
270  - (this->alpha_*this->rho_*this->nu())
271  *dev(twoSymm(fvc::grad(this->U_)))
272  );
273 }
274 
275 
276 template<class BasicMomentumTransportModel>
279 (
280  volVectorField& U
281 ) const
282 {
283  return
284  (
285  fvc::div
286  (
287  this->alpha_*this->rho_*this->nuM_*fvc::grad(U)
288  )
289  + fvc::div(this->alpha_*this->rho_*sigma_)
290  - fvc::div(this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(U))))
291  - fvm::laplacian(this->alpha_*this->rho_*nu0(), U)
292  );
293 }
294 
295 
296 template<class BasicMomentumTransportModel>
299 (
300  const volScalarField& rho,
301  volVectorField& U
302 ) const
303 {
304  return
305  (
306  fvc::div
307  (
308  this->alpha_*rho*this->nuM_*fvc::grad(U)
309  )
310  + fvc::div(this->alpha_*rho*sigma_)
311  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
312  - fvm::laplacian(this->alpha_*rho*nu0(), U)
313  );
314 }
315 
316 
317 template<class BasicMomentumTransportModel>
319 {
320  // Local references
321  const alphaField& alpha = this->alpha_;
322  const rhoField& rho = this->rho_;
323  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
324  const volVectorField& U = this->U_;
325  fv::options& fvOptions(fv::options::New(this->mesh_));
326 
328 
329  tmp<volTensorField> tgradU(fvc::grad(U));
330  const volTensorField& gradU = tgradU();
331 
332  forAll(lambdas_, modei)
333  {
334  volSymmTensorField& sigma = nModes_ == 1 ? sigma_ : sigmas_[modei];
335 
337  (
338  IOobject
339  (
341  (
342  "rLambda"
343  + (nModes_ == 1 ? word::null : name(modei)),
344  this->alphaRhoPhi_.group()
345  ),
346  this->runTime_.constant(),
347  this->mesh_
348  ),
349  1/lambdas_[modei]
350  );
351 
352  // Note sigma is positive on lhs of momentum eqn
353  const volSymmTensorField P
354  (
355  twoSymm(sigma & gradU)
356  - nuM_*rLambda*twoSymm(gradU)
357  );
358 
359  // Viscoelastic stress equation
360  fvSymmTensorMatrix sigmaEqn
361  (
362  fvm::ddt(alpha, rho, sigma)
363  + fvm::div
364  (
365  alphaRhoPhi,
366  sigma,
367  "div(" + alphaRhoPhi.name() + ',' + sigma_.name() + ')'
368  )
369  ==
370  alpha*rho*P
371  + sigmaSource(modei, sigma)
372  + fvOptions(alpha, rho, sigma)
373  );
374 
375  sigmaEqn.relax();
376  fvOptions.constrain(sigmaEqn);
377  sigmaEqn.solve("sigma");
378  fvOptions.correct(sigma);
379  }
380 
381  if (sigmas_.size())
382  {
383  volSymmTensorField sigmaSum("sigmaSum", sigmas_[0]);
384 
385  for (label modei = 1; modei<sigmas_.size(); modei++)
386  {
387  sigmaSum += sigmas_[modei];
388  }
389 
390  sigma_ == sigmaSum;
391  }
392 }
393 
394 
395 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
396 
397 } // End namespace laminarModels
398 } // End namespace Foam
399 
400 // ************************************************************************* //
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
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:303
fv::options & fvOptions
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
Maxwell(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: Maxwell.C:105
PtrList< dimensionedScalar > readModeCoefficients(const word &name, const dimensionSet &dims) const
Definition: Maxwell.C:42
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< symmTensor >> &)
Return a temporary field constructed from name,.
const dimensionSet dimViscosity
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
Definition: Maxwell.C:279
BasicMomentumTransportModel::rhoField rhoField
Definition: laminarModel.H:77
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
BasicMomentumTransportModel::alphaField alphaField
Definition: laminarModel.H:76
Finite-volume options.
Definition: fvOptions.H:52
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
Definition: Maxwell.C:264
BasicMomentumTransportModel::transportModel transportModel
Definition: laminarModel.H:78
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
Dimension set for the base types.
Definition: dimensionSet.H:120
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
virtual tmp< fvSymmTensorMatrix > sigmaSource(const label modei, volSymmTensorField &sigma) const
Definition: Maxwell.C:92
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
static const word null
An empty word.
Definition: word.H:77
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
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:232
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:521
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: Maxwell.C:318
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
Definition: Maxwell.C:256
U
Definition: pEqn.H:72
virtual void correct()
Correct the laminar transport.
Definition: laminarModel.C:271
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
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
messageStream Info
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
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
Namespace for OpenFOAM.