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-2024 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 "fvModels.H"
28 #include "fvConstraints.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace laminarModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicMomentumTransportModel>
41 PtrList<dimensionedScalar>
43 (
44  const word& name,
45  const dimensionSet& dims
46 ) const
47 {
48  PtrList<dimensionedScalar> modeCoeffs(nModes_);
49 
50  if (modeCoefficients_.size())
51  {
52  if (this->coeffDict().found(name))
53  {
54  IOWarningInFunction(this->coeffDict())
55  << "Using 'modes' list, '" << name << "' entry will be ignored."
56  << endl;
57  }
58 
59  forAll(modeCoefficients_, modei)
60  {
61  modeCoeffs.set
62  (
63  modei,
65  (
66  name,
67  dims,
68  modeCoefficients_[modei].lookup(name)
69  )
70  );
71  }
72  }
73  else
74  {
75  modeCoeffs.set
76  (
77  0,
79  (
80  name,
81  dims,
82  this->coeffDict().lookup(name)
83  )
84  );
85  }
86 
87  return modeCoeffs;
88 }
89 
90 
91 template<class BasicMomentumTransportModel>
93 (
94  const label modei,
96 ) const
97 {
98  return -fvm::Sp(this->alpha_*this->rho_/lambdas_[modei], sigma);
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103 
104 template<class BasicMomentumTransportModel>
106 (
107  const alphaField& alpha,
108  const rhoField& rho,
109  const volVectorField& U,
110  const surfaceScalarField& alphaRhoPhi,
111  const surfaceScalarField& phi,
112  const viscosity& viscosity,
113  const word& type
114 )
115 :
116  laminarModel<BasicMomentumTransportModel>
117  (
118  type,
119  alpha,
120  rho,
121  U,
122  alphaRhoPhi,
123  phi,
124  viscosity
125  ),
126 
127  modeCoefficients_
128  (
129  this->coeffDict().found("modes")
131  (
132  this->coeffDict().lookup("modes")
133  )
134  : PtrList<dictionary>()
135  ),
136 
137  nModes_(modeCoefficients_.size() ? modeCoefficients_.size() : 1),
138 
139  nuM_("nuM", dimKinematicViscosity, this->coeffDict().lookup("nuM")),
140 
141  lambdas_(readModeCoefficients("lambda", dimTime)),
142 
143  sigma_
144  (
145  IOobject
146  (
147  this->groupName("sigma"),
148  this->runTime_.name(),
149  this->mesh_,
150  IOobject::MUST_READ,
151  IOobject::AUTO_WRITE
152  ),
153  this->mesh_
154  )
155 {
156  if (nModes_ > 1)
157  {
158  sigmas_.setSize(nModes_);
159 
160  forAll(sigmas_, modei)
161  {
163  (
164  this->groupName("sigma" + name(modei)),
165  this->runTime_.name(),
166  this->mesh_,
168  );
169 
170  // Check if mode field exists and can be read
171  if (header.headerOk())
172  {
173  Info<< " Reading mode stress field "
174  << header.name() << endl;
175 
176  sigmas_.set
177  (
178  modei,
180  (
181  IOobject
182  (
183  header.name(),
184  this->runTime_.name(),
185  this->mesh_,
188  ),
189  this->mesh_
190  )
191  );
192  }
193  else
194  {
195  sigmas_.set
196  (
197  modei,
199  (
200  IOobject
201  (
202  header.name(),
203  this->runTime_.name(),
204  this->mesh_,
207  ),
208  sigma_
209  )
210  );
211  }
212  }
213  }
214 }
215 
216 
217 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218 
219 template<class BasicMomentumTransportModel>
221 {
223  {
224  if (modeCoefficients_.size())
225  {
226  this->coeffDict().lookup("modes") >> modeCoefficients_;
227  }
228 
229  nuM_.read(this->coeffDict());
230 
231  lambdas_ = readModeCoefficients("lambda", dimTime);
232 
233  return true;
234  }
235  else
236  {
237  return false;
238  }
239 }
240 
241 
242 template<class BasicMomentumTransportModel>
244 {
245  return volScalarField::New
246  (
247  this->groupName("nuEff"),
248  this->nu()
249  );
250 }
251 
252 
253 template<class BasicMomentumTransportModel>
255 (
256  const label patchi
257 ) const
258 {
259  return this->nu(patchi);
260 }
261 
262 
263 template<class BasicMomentumTransportModel>
265 {
266  const surfaceVectorField nf(this->mesh().nf());
267 
269  (
270  this->groupName("devTau"),
272  (
273  nf,
274  this->alpha_*this->rho_*this->nuM_*fvc::grad(this->U_)
275  )
277  (
278  nf,
279  this->alpha_*this->rho_*sigma_
280  )
282  (
283  nf,
284  this->alpha_*this->rho_*this->nu()*dev2(T(fvc::grad(this->U_)))
285  )
286  - fvc::interpolate(this->alpha_*this->rho_*nu0())*fvc::snGrad(this->U_)
287  );
288 }
289 
290 
291 template<class BasicMomentumTransportModel>
292 template<class RhoFieldType>
295 (
296  const RhoFieldType& rho,
298 ) const
299 {
300  return
301  (
302  fvc::div(this->alpha_*rho*this->nuM_*fvc::grad(U))
303  + fvc::div(this->alpha_*rho*sigma_)
304  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
305  - fvm::laplacian(this->alpha_*rho*nu0(), U)
306  );
307 }
308 
309 
310 template<class BasicMomentumTransportModel>
312 (
314 ) const
315 {
316  return DivDevTau(this->rho_, U);
317 }
318 
319 
320 template<class BasicMomentumTransportModel>
323 (
324  const volScalarField& rho,
326 ) const
327 {
328  return DivDevTau(rho, U);
329 }
330 
331 
332 template<class BasicMomentumTransportModel>
334 {
335  // Local references
336  const alphaField& alpha = this->alpha_;
337  const rhoField& rho = this->rho_;
338  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
339  const volVectorField& U = this->U_;
340  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
342  (
343  Foam::fvConstraints::New(this->mesh_)
344  );
345 
347 
349  const volTensorField& gradU = tgradU();
350 
351  forAll(lambdas_, modei)
352  {
353  volSymmTensorField& sigma = nModes_ == 1 ? sigma_ : sigmas_[modei];
354 
356  (
357  IOobject
358  (
359  this->groupName
360  (
361  "rLambda"
362  + (nModes_ == 1 ? word::null : name(modei))
363  ),
364  this->runTime_.constant(),
365  this->mesh_
366  ),
367  1/lambdas_[modei]
368  );
369 
370  // Note sigma is positive on lhs of momentum eqn
371  const volSymmTensorField P
372  (
373  twoSymm(sigma & gradU)
374  - nuM_*rLambda*twoSymm(gradU)
375  );
376 
377  // Viscoelastic stress equation
378  fvSymmTensorMatrix sigmaEqn
379  (
381  + fvm::div
382  (
383  alphaRhoPhi,
384  sigma,
385  "div(" + alphaRhoPhi.name() + ',' + sigma_.name() + ')'
386  )
387  ==
388  alpha*rho*P
389  + sigmaSource(modei, sigma)
391  );
392 
393  sigmaEqn.relax();
394  fvConstraints.constrain(sigmaEqn);
395  sigmaEqn.solve("sigma");
397  }
398 
399  if (sigmas_.size())
400  {
401  volSymmTensorField sigmaSum("sigmaSum", sigmas_[0]);
402 
403  for (label modei = 1; modei<sigmas_.size(); modei++)
404  {
405  sigmaSum += sigmas_[modei];
406  }
407 
408  sigma_ == sigmaSum;
409  }
410 }
411 
412 
413 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
414 
415 } // End namespace laminarModels
416 } // End namespace Foam
417 
418 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:307
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Dimension set for the base types.
Definition: dimensionSet.H:125
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:603
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
Templated abstract base class for laminar transport models.
Definition: laminarModel.H:52
BasicMomentumTransportModel::alphaField alphaField
Definition: laminarModel.H:67
virtual void correct()
Predict the laminar viscosity.
Definition: laminarModel.C:274
BasicMomentumTransportModel::rhoField rhoField
Definition: laminarModel.H:68
Generalised Maxwell model for viscoelasticity using the upper-convected time derivative of the stress...
Definition: Maxwell.H:79
volSymmTensorField sigma_
Single or mode sum viscoelastic stress.
Definition: Maxwell.H:109
virtual tmp< fvSymmTensorMatrix > sigmaSource(const label modei, volSymmTensorField &sigma) const
Definition: Maxwell.C:93
virtual void correct()
Correct the Maxwell stress.
Definition: Maxwell.C:333
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity, i.e. the laminar viscosity.
Definition: Maxwell.C:243
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
Definition: Maxwell.C:312
virtual tmp< surfaceVectorField > devTau() const
Return the effective surface stress.
Definition: Maxwell.C:264
PtrList< volSymmTensorField > sigmas_
Mode viscoelastic stresses.
Definition: Maxwell.H:112
Maxwell(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: Maxwell.C:106
PtrList< dimensionedScalar > readModeCoefficients(const word &name, const dimensionSet &dims) const
Definition: Maxwell.C:43
virtual bool read()
Read model coefficients if they have changed.
Definition: Maxwell.C:220
A class for managing temporary objects.
Definition: tmp.H:55
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:530
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate(const surfaceVectorField &Sf, const VolField< Type > &tvf)
Interpolate field onto faces.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
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
void dev2(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimKinematicViscosity
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
messageStream Info
const dimensionSet dimTime
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488