kineticTheoryModel.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) 2011-2023 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 "kineticTheoryModel.H"
27 #include "mathematicalConstants.H"
28 #include "phaseSystem.H"
29 #include "fvcDdt.H"
30 #include "fvcSup.H"
31 #include "fvModels.H"
32 #include "fvConstraints.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 const Foam::phaseModel&
37 Foam::RASModels::kineticTheoryModel::continuousPhase() const
38 {
39  const phaseSystem& fluid = phase_.fluid();
40 
41  if (continuousPhaseName_ == word::null)
42  {
43  if (fluid.movingPhases().size() != 2)
44  {
45  FatalIOErrorInFunction(coeffDict_)
46  << "Continuous phase name must be specified "
47  << "when there are more than two moving phases."
48  << exit(FatalIOError);
49  }
50 
51  forAll(fluid.movingPhases(), movingPhasei)
52  {
53  const phaseModel& otherPhase = fluid.movingPhases()[movingPhasei];
54 
55  if (&otherPhase != &phase_)
56  {
57  return otherPhase;
58  }
59  }
60  }
61 
62  return fluid.phases()[continuousPhaseName_];
63 };
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
69 (
70  const volScalarField& alpha,
71  const volScalarField& rho,
72  const volVectorField& U,
73  const surfaceScalarField& alphaRhoPhi,
74  const surfaceScalarField& phi,
75  const viscosity& viscosity,
76  const word& type
77 )
78 :
79  eddyViscosity<RASModel<phaseCompressible::momentumTransportModel>>
80  (
81  type,
82  alpha,
83  rho,
84  U,
85  alphaRhoPhi,
86  phi,
87  viscosity
88  ),
89 
90  phase_(refCast<const phaseModel>(viscosity)),
91 
92  continuousPhaseName_
93  (
94  coeffDict_.lookupOrDefault("continuousPhase", word::null)
95  ),
96 
97  viscosityModel_
98  (
99  kineticTheoryModels::viscosityModel::New
100  (
101  coeffDict_
102  )
103  ),
104  conductivityModel_
105  (
106  kineticTheoryModels::conductivityModel::New
107  (
108  coeffDict_
109  )
110  ),
111  radialModel_
112  (
113  kineticTheoryModels::radialModel::New
114  (
115  coeffDict_
116  )
117  ),
118  granularPressureModel_
119  (
120  kineticTheoryModels::granularPressureModel::New
121  (
122  coeffDict_
123  )
124  ),
125  frictionalStressModel_
126  (
127  kineticTheoryModels::frictionalStressModel::New
128  (
129  coeffDict_
130  )
131  ),
132 
133  equilibrium_(coeffDict_.lookup("equilibrium")),
134  e_("e", dimless, coeffDict_),
135  alphaMinFriction_
136  (
137  "alphaMinFriction",
138  dimless,
139  coeffDict_
140  ),
141  residualAlpha_
142  (
143  "residualAlpha",
144  dimless,
145  coeffDict_
146  ),
147 
148  maxNut_
149  (
150  "maxNut",
151  dimensionSet(0, 2, -1, 0, 0),
152  coeffDict_.lookupOrDefault<scalar>("maxNut", 1000)
153  ),
154 
155  Theta_
156  (
157  IOobject
158  (
159  IOobject::groupName("Theta", phase_.name()),
160  U.time().name(),
161  U.mesh(),
162  IOobject::MUST_READ,
163  IOobject::AUTO_WRITE
164  ),
165  U.mesh()
166  ),
167 
168  lambda_
169  (
170  IOobject
171  (
172  IOobject::groupName(typedName("lambda"), phase_.name()),
173  U.time().name(),
174  U.mesh(),
175  IOobject::NO_READ,
176  IOobject::NO_WRITE
177  ),
178  U.mesh(),
179  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), 0)
180  ),
181 
182  gs0_
183  (
184  IOobject
185  (
186  IOobject::groupName(typedName("gs0"), phase_.name()),
187  U.time().name(),
188  U.mesh(),
189  IOobject::NO_READ,
190  IOobject::NO_WRITE
191  ),
192  U.mesh(),
193  dimensionedScalar(dimensionSet(0, 0, 0, 0, 0), 0)
194  ),
195 
196  kappa_
197  (
198  IOobject
199  (
200  IOobject::groupName(typedName("kappa"), phase_.name()),
201  U.time().name(),
202  U.mesh(),
203  IOobject::NO_READ,
204  IOobject::NO_WRITE
205  ),
206  U.mesh(),
207  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), 0)
208  ),
209 
210  nuFric_
211  (
212  IOobject
213  (
214  IOobject::groupName(typedName("nuFric"), phase_.name()),
215  U.time().name(),
216  U.mesh(),
217  IOobject::NO_READ,
218  IOobject::AUTO_WRITE
219  ),
220  U.mesh(),
221  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), 0)
222  )
223 {
224  if (type == typeName)
225  {
226  printCoeffs(type);
227  }
228 }
229 
230 
231 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
232 
234 {}
235 
236 
237 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
238 
240 {
241  if
242  (
244  read()
245  )
246  {
247  coeffDict().lookup("equilibrium") >> equilibrium_;
248  e_.readIfPresent(coeffDict());
249  alphaMinFriction_.readIfPresent(coeffDict());
250 
251  viscosityModel_->read();
252  conductivityModel_->read();
253  radialModel_->read();
254  granularPressureModel_->read();
255  frictionalStressModel_->read();
256 
257  return true;
258  }
259  else
260  {
261  return false;
262  }
263 }
264 
265 
268 {
270  return nut_;
271 }
272 
273 
276 {
278  return nut_;
279 }
280 
281 
284 {
286  return nut_;
287 }
288 
289 
292 {
294  (
296  (
297  IOobject::groupName("R", U_.group()),
298  - (nut_ + nuFric_)*dev(twoSymm(fvc::grad(U_)))
299  - (lambda_*fvc::div(phi_))*symmTensor::I
300  )
301  );
302 }
303 
304 
307 {
308  const volScalarField& rho = phase_.rho();
309 
310  tmp<volScalarField> tpPrime
311  (
313  (
314  IOobject::groupName("pPrime", U_.group()),
315  Theta_
316  *granularPressureModel_->granularPressureCoeffPrime
317  (
318  alpha_,
319  radialModel_->g0
320  (
321  alpha_,
322  alphaMinFriction_,
323  phase_.alphaMax()
324  ),
325  radialModel_->g0prime
326  (
327  alpha_,
328  alphaMinFriction_,
329  phase_.alphaMax()
330  ),
331  rho,
332  e_
333  )
334  + frictionalStressModel_->frictionalPressurePrime
335  (
336  phase_,
337  alphaMinFriction_,
338  phase_.alphaMax()
339  )
340  )
341  );
342 
343  volScalarField::Boundary& bpPrime = tpPrime.ref().boundaryFieldRef();
344 
345  forAll(bpPrime, patchi)
346  {
347  if (!bpPrime[patchi].coupled())
348  {
349  bpPrime[patchi] == 0;
350  }
351  }
352 
353  return tpPrime;
354 }
355 
356 
359 {
361  (
362  IOobject::groupName("pPrimef", U_.group()),
363  fvc::interpolate(pPrime())
364  );
365 }
366 
367 
370 {
372  (
374  (
375  IOobject::groupName("devTau", U_.group()),
376  - (rho_*(nut_ + nuFric_))
377  *dev(twoSymm(fvc::grad(U_)))
378  - ((rho_*lambda_)*fvc::div(phi_))*symmTensor::I
379  )
380  );
381 }
382 
383 
386 (
388 ) const
389 {
390  return
391  (
392  - fvm::laplacian(rho_*(nut_ + nuFric_), U)
393  - fvc::div
394  (
395  (rho_*(nut_ + nuFric_))*dev2(T(fvc::grad(U)))
396  + ((rho_*lambda_)*fvc::div(phi_))
398  "divDevTau(" + U_.name() + ')'
399  )
400  );
401 }
402 
403 
405 {
406  // Local references
407  const volScalarField alpha(max(alpha_, scalar(0)));
408  const phaseSystem& fluid = phase_.fluid();
409  const phaseModel& continuousPhase = this->continuousPhase();
410  const volScalarField& rho = phase_.rho();
411  const surfaceScalarField& alphaRhoPhi = alphaRhoPhi_;
412  const volVectorField& U = U_;
413 
414  tmp<volVectorField> tUc(continuousPhase.U());
415  const volVectorField& Uc = tUc();
416 
417  const scalar sqrtPi = sqrt(constant::mathematical::pi);
418  const dimensionedScalar ThetaSmall("ThetaSmall", Theta_.dimensions(), 1e-6);
419  const dimensionedScalar ThetaSmallSqrt(sqrt(ThetaSmall));
420 
421  tmp<volScalarField> tda(phase_.d());
422  const volScalarField& da = tda();
423 
424  tmp<volTensorField> tgradU(fvc::grad(U_));
425  const volTensorField& gradU(tgradU());
426  const volSymmTensorField D(symm(gradU));
427 
428  // Calculating the radial distribution function
429  gs0_ = radialModel_->g0(alpha, alphaMinFriction_, phase_.alphaMax());
430 
431  if (!equilibrium_)
432  {
433  // Particle viscosity (Table 3.2, p.47)
434  nut_ = viscosityModel_->nu(alpha, Theta_, gs0_, rho, da, e_);
435 
436  const volScalarField ThetaSqrt("sqrtTheta", sqrt(Theta_));
437 
438  // Bulk viscosity p. 45 (Lun et al. 1984).
439  lambda_ = (4.0/3.0)*sqr(alpha)*da*gs0_*(1 + e_)*ThetaSqrt/sqrtPi;
440 
441  // Stress tensor, Definitions, Table 3.1, p. 43
442  const volSymmTensorField tau
443  (
444  rho*(2*nut_*D + (lambda_ - (2.0/3.0)*nut_)*tr(D)*I)
445  );
446 
447  // Dissipation (Eq. 3.24, p.50)
448  const volScalarField gammaCoeff
449  (
450  "gammaCoeff",
451  12*(1 - sqr(e_))
452  *max(sqr(alpha), residualAlpha_)
453  *rho*gs0_*(1.0/da)*ThetaSqrt/sqrtPi
454  );
455 
456  // Drag
457  const dispersedPhaseInterface interface(phase_, continuousPhase);
458  const volScalarField beta
459  (
460  fluid.foundInterfacialModel<dragModel>(interface)
461  ? fluid.lookupInterfacialModel<dragModel>(interface).K()
463  (
464  "beta",
465  phase_.mesh(),
467  )
468  );
469 
470  // Eq. 3.25, p. 50 Js = J1 - J2
471  const volScalarField J1("J1", 3*beta);
472  const volScalarField J2
473  (
474  "J2",
475  0.25*sqr(beta)*da*magSqr(U - Uc)
476  /(
477  max(alpha, residualAlpha_)*rho
478  *sqrtPi*(ThetaSqrt + ThetaSmallSqrt)
479  )
480  );
481 
482  // particle pressure - coefficient in front of Theta (Eq. 3.22, p. 45)
483  const volScalarField PsCoeff
484  (
485  granularPressureModel_->granularPressureCoeff
486  (
487  alpha,
488  gs0_,
489  rho,
490  e_
491  )
492  );
493 
494  // 'thermal' conductivity (Table 3.3, p. 49)
495  kappa_ = conductivityModel_->kappa(alpha, Theta_, gs0_, rho, da, e_);
496 
499  (
501  );
502 
503  // Construct the granular temperature equation (Eq. 3.20, p. 44)
504  // NB. note that there are two typos in Eq. 3.20:
505  // Ps should be without grad
506  // the laplacian has the wrong sign
507  fvScalarMatrix ThetaEqn
508  (
509  1.5*
510  (
511  fvm::ddt(alpha, rho, Theta_)
512  + fvm::div(alphaRhoPhi, Theta_)
513  - fvc::Sp(fvc::ddt(alpha, rho) + fvc::div(alphaRhoPhi), Theta_)
514  )
515  - fvm::laplacian(kappa_, Theta_, "laplacian(kappa,Theta)")
516  ==
517  - fvm::SuSp((PsCoeff*I) && gradU, Theta_)
518  + (tau && gradU)
519  + fvm::Sp(-gammaCoeff, Theta_)
520  + fvm::Sp(-J1, Theta_)
521  + fvm::Sp(J2/(Theta_ + ThetaSmall), Theta_)
522  + fvModels.source(alpha, rho, Theta_)
523  );
524 
525  ThetaEqn.relax();
526  fvConstraints.constrain(ThetaEqn);
527  ThetaEqn.solve();
528  fvConstraints.constrain(Theta_);
529  }
530  else
531  {
532  // Equilibrium => dissipation == production
533  // Eq. 4.14, p.82
534  const volScalarField K1("K1", 2*(1 + e_)*rho*gs0_);
535  const volScalarField K3
536  (
537  "K3",
538  0.5*da*rho*
539  (
540  (sqrtPi/(3*(3.0 - e_)))
541  *(1 + 0.4*(1 + e_)*(3*e_ - 1)*alpha*gs0_)
542  +1.6*alpha*gs0_*(1 + e_)/sqrtPi
543  )
544  );
545 
546  const volScalarField K2
547  (
548  "K2",
549  4*da*rho*(1 + e_)*alpha*gs0_/(3*sqrtPi) - 2*K3/3.0
550  );
551 
552  const volScalarField K4("K4", 12*(1 - sqr(e_))*rho*gs0_/(da*sqrtPi));
553 
554  const volScalarField trD
555  (
556  "trD",
557  alpha/(alpha + residualAlpha_)
558  *fvc::div(phi_)
559  );
560  const volScalarField tr2D("tr2D", sqr(trD));
561  const volScalarField trD2("trD2", tr(D & D));
562 
563  const volScalarField t1("t1", K1*alpha + rho);
564  const volScalarField l1("l1", -t1*trD);
565  const volScalarField l2("l2", sqr(t1)*tr2D);
566  const volScalarField l3
567  (
568  "l3",
569  4.0
570  *K4
571  *alpha
572  *(2*K3*trD2 + K2*tr2D)
573  );
574 
575  Theta_ = sqr
576  (
577  (l1 + sqrt(l2 + l3))
578  /(2*max(alpha, residualAlpha_)*K4)
579  );
580 
581  kappa_ = conductivityModel_->kappa(alpha, Theta_, gs0_, rho, da, e_);
582  }
583 
584  Theta_.max(0);
585  Theta_.min(100);
586 
587  {
588  // particle viscosity (Table 3.2, p.47)
589  nut_ = viscosityModel_->nu(alpha, Theta_, gs0_, rho, da, e_);
590 
591  const volScalarField ThetaSqrt("sqrtTheta", sqrt(Theta_));
592 
593  // Bulk viscosity p. 45 (Lun et al. 1984).
594  lambda_ = (4.0/3.0)*sqr(alpha)*da*gs0_*(1 + e_)*ThetaSqrt/sqrtPi;
595 
596  // Frictional pressure
597  const volScalarField pf
598  (
599  frictionalStressModel_->frictionalPressure
600  (
601  phase_,
602  alphaMinFriction_,
603  phase_.alphaMax()
604  )
605  );
606 
607  // Limit viscosity
608  nut_.min(maxNut_);
609 
610  nuFric_ = min
611  (
612  frictionalStressModel_->nu
613  (
614  phase_,
615  alphaMinFriction_,
616  phase_.alphaMax(),
617  pf/rho,
618  D
619  ),
620  maxNut_ - nut_
621  );
622  }
623 
624  if (debug)
625  {
626  Info<< typeName << ':' << nl
627  << " max(Theta) = " << max(Theta_).value() << nl
628  << " max(nut) = " << max(nut_).value() << endl;
629  }
630 }
631 
632 
633 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
#define K3
Definition: SHA1.C:169
#define K1
Definition: SHA1.C:167
#define K4
Definition: SHA1.C:170
#define K2
Definition: SHA1.C:168
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricBoundaryField class.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
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
static word groupName(Name name, const word &group)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< volSymmTensorField > devTau() const
Return the effective stress tensor.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volSymmTensorField > sigma() const
Return the stress tensor [m^2/s^2].
kineticTheoryModel(const volScalarField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual bool read()
Re-read model coefficients if they have changed.
static const SymmTensor I
Definition: SymmTensor.H:72
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Dimension set for the base types.
Definition: dimensionSet.H:122
const word & name() const
Return const reference to name.
Class to represent a interface between phases where one phase is considered dispersed within the othe...
Model for drag between phases.
Definition: dragModel.H:55
virtual tmp< volScalarField > K() const =0
Return the drag coefficient K.
static const dimensionSet dimK
Coefficient dimensions.
Definition: dragModel.H:81
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite volume constraints.
Definition: fvConstraints.H:62
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:604
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
Finite volume models.
Definition: fvModels.H:65
Abstract base class for turbulence models (RAS, LES and laminar).
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:127
virtual tmp< volVectorField > U() const =0
Return the velocity.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
bool foundInterfacialModel(const phaseInterface &interface) const
Check availability of a sub model for a given interface.
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
An abstract base class for Newtonian viscosity models.
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))
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
Calculate the first temporal derivative.
Calculate the field for explicit evaluation of implicit and explicit sources.
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))
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
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< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
const doubleScalar e
Definition: doubleScalar.H:105
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor sqr(const dimensionedVector &dv)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:111
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimless
messageStream Info
static const Identity< scalar > I
Definition: Identity.H:93
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:149
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488