ThermalPhaseChangePhaseSystem.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2019 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "fvcVolumeIntegrate.H"
31 #include "fvmSup.H"
32 
33 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
34 
35 template<class BasePhaseSystem>
38 (
39  const phasePairKey& key
40 ) const
41 {
42  if (!iDmdt_.found(key))
43  {
44  return phaseSystem::dmdt(key);
45  }
46 
47  const scalar dmdtSign(Pair<word>::compare(iDmdt_.find(key).key(), key));
48 
49  return dmdtSign**iDmdt_[key];
50 }
51 
52 
53 template<class BasePhaseSystem>
56 (
57  const phasePairKey& key
58 ) const
59 {
60  if (!wDmdt_.found(key))
61  {
62  return phaseSystem::dmdt(key);
63  }
64 
65  const scalar dmdtSign(Pair<word>::compare(wDmdt_.find(key).key(), key));
66 
67  return dmdtSign**wDmdt_[key];
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72 
73 template<class BasePhaseSystem>
76 (
77  const fvMesh& mesh
78 )
79 :
80  BasePhaseSystem(mesh),
81  volatile_(this->template lookupOrDefault<word>("volatile", "none")),
82  saturationModel_
83  (
84  saturationModel::New(this->subDict("saturationModel"), mesh)
85  ),
86  phaseChange_(this->lookup("phaseChange"))
87 {
88 
90  (
92  this->phasePairs_,
93  phasePairIter
94  )
95  {
96  const phasePair& pair(phasePairIter());
97 
98  if (pair.ordered())
99  {
100  continue;
101  }
102 
103  // Initially assume no mass transfer
104  iDmdt_.set
105  (
106  pair,
107  new volScalarField
108  (
109  IOobject
110  (
111  IOobject::groupName("iDmdt", pair.name()),
112  this->mesh().time().timeName(),
113  this->mesh(),
114  IOobject::READ_IF_PRESENT,
115  IOobject::AUTO_WRITE
116  ),
117  this->mesh(),
119  )
120  );
121 
122  // Initially assume no mass transfer
123  wDmdt_.set
124  (
125  pair,
126  new volScalarField
127  (
128  IOobject
129  (
130  IOobject::groupName("wDmdt", pair.name()),
131  this->mesh().time().timeName(),
132  this->mesh(),
133  IOobject::READ_IF_PRESENT,
134  IOobject::AUTO_WRITE
135  ),
136  this->mesh(),
138  )
139  );
140 
141  // Initially assume no mass transfer
142  wMDotL_.set
143  (
144  pair,
145  new volScalarField
146  (
147  IOobject
148  (
149  IOobject::groupName("wMDotL", pair.name()),
150  this->mesh().time().timeName(),
151  this->mesh(),
152  IOobject::READ_IF_PRESENT,
153  IOobject::AUTO_WRITE
154  ),
155  this->mesh(),
157  )
158  );
159  }
160 }
161 
162 
163 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
164 
165 template<class BasePhaseSystem>
168 {}
169 
170 
171 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
172 
173 template<class BasePhaseSystem>
176 {
177  return saturationModel_();
178 }
179 
180 
181 template<class BasePhaseSystem>
184 (
185  const phasePairKey& key
186 ) const
187 {
188  return BasePhaseSystem::dmdt(key) + this->iDmdt(key) + this->wDmdt(key);
189 }
190 
191 
192 template<class BasePhaseSystem>
195 {
196  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
197 
198  forAllConstIter(iDmdtTable, iDmdt_, iDmdtIter)
199  {
200  const phasePair& pair = this->phasePairs_[iDmdtIter.key()];
201  const volScalarField& iDmdt = *iDmdtIter();
202 
203  this->addField(pair.phase1(), "dmdt", iDmdt, dmdts);
204  this->addField(pair.phase2(), "dmdt", - iDmdt, dmdts);
205  }
206 
207  forAllConstIter(wDmdtTable, wDmdt_, wDmdtIter)
208  {
209  const phasePair& pair = this->phasePairs_[wDmdtIter.key()];
210  const volScalarField& wDmdt = *wDmdtIter();
211 
212  this->addField(pair.phase1(), "dmdt", wDmdt, dmdts);
213  this->addField(pair.phase2(), "dmdt", - wDmdt, dmdts);
214  }
215 
216  return dmdts;
217 }
218 
219 
220 template<class BasePhaseSystem>
223 {
225  BasePhaseSystem::heatTransfer();
226 
227  phaseSystem::heatTransferTable& eqns = eqnsPtr();
228 
229  // Add boundary term
231  (
233  this->phasePairs_,
234  phasePairIter
235  )
236  {
237  if (this->wMDotL_.found(phasePairIter.key()))
238  {
239  const phasePair& pair(phasePairIter());
240 
241  if (pair.ordered())
242  {
243  continue;
244  }
245 
246  const phaseModel& phase1 = pair.phase1();
247  const phaseModel& phase2 = pair.phase2();
248 
249  *eqns[phase1.name()] += negPart(*this->wMDotL_[pair]);
250  *eqns[phase2.name()] -= posPart(*this->wMDotL_[pair]);
251 
252  if
253  (
254  phase1.thermo().he().member() == "e"
255  || phase2.thermo().he().member() == "e"
256  )
257  {
258  const volScalarField dmdt
259  (
260  this->iDmdt(pair) + this->wDmdt(pair)
261  );
262 
263  if (phase1.thermo().he().member() == "e")
264  {
265  *eqns[phase1.name()] +=
266  phase1.thermo().p()*dmdt/phase1.thermo().rho();
267  }
268 
269  if (phase2.thermo().he().member() == "e")
270  {
271  *eqns[phase2.name()] -=
272  phase2.thermo().p()*dmdt/phase2.thermo().rho();
273  }
274  }
275  }
276  }
277 
278  return eqnsPtr;
279 }
280 
281 
282 template<class BasePhaseSystem>
285 {
288 
289  phaseSystem::massTransferTable& eqns = eqnsPtr();
290 
292  (
294  this->phasePairs_,
295  phasePairIter
296  )
297  {
298  const phasePair& pair(phasePairIter());
299 
300  if (pair.ordered())
301  {
302  continue;
303  }
304 
305  const phaseModel& phase = pair.phase1();
306  const phaseModel& otherPhase = pair.phase2();
307 
308  const PtrList<volScalarField>& Yi = phase.Y();
309 
310  forAll(Yi, i)
311  {
312  if (Yi[i].member() != volatile_)
313  {
314  continue;
315  }
316 
317  const word name
318  (
319  IOobject::groupName(volatile_, phase.name())
320  );
321 
322  const word otherName
323  (
324  IOobject::groupName(volatile_, otherPhase.name())
325  );
326 
327  // Note that the phase YiEqn does not contain a continuity error
328  // term, so these additions represent the entire mass transfer
329 
330  const volScalarField dmdt(this->iDmdt(pair) + this->wDmdt(pair));
331 
332  *eqns[name] += dmdt;
333  *eqns[otherName] -= dmdt;
334  }
335  }
336 
337  return eqnsPtr;
338 }
339 
340 
341 template<class BasePhaseSystem>
342 void
344 {
346  alphatPhaseChangeWallFunction;
347 
349  (
350  typename BasePhaseSystem::heatTransferModelTable,
351  this->heatTransferModels_,
352  heatTransferModelIter
353  )
354  {
355  const phasePair& pair
356  (
357  this->phasePairs_[heatTransferModelIter.key()]
358  );
359 
360  const phaseModel& phase1 = pair.phase1();
361  const phaseModel& phase2 = pair.phase2();
362 
363  const volScalarField& T1(phase1.thermo().T());
364  const volScalarField& T2(phase2.thermo().T());
365 
366  const volScalarField& he1(phase1.thermo().he());
367  const volScalarField& he2(phase2.thermo().he());
368 
369  const volScalarField& p(phase1.thermo().p());
370 
371  volScalarField& iDmdt(*this->iDmdt_[pair]);
372  volScalarField& Tf(*this->Tf_[pair]);
373 
374  const volScalarField Tsat(saturationModel_->Tsat(phase1.thermo().p()));
375 
376  volScalarField hf1
377  (
378  he1.member() == "e"
379  ? phase1.thermo().he(p, Tsat) + p/phase1.rho()
380  : phase1.thermo().he(p, Tsat)
381  );
382  volScalarField hf2
383  (
384  he2.member() == "e"
385  ? phase2.thermo().he(p, Tsat) + p/phase2.rho()
386  : phase2.thermo().he(p, Tsat)
387  );
388 
389  volScalarField h1
390  (
391  he1.member() == "e"
392  ? he1 + p/phase1.rho()
393  : tmp<volScalarField>(he1)
394  );
395 
396  volScalarField h2
397  (
398  he2.member() == "e"
399  ? he2 + p/phase2.rho()
401  );
402 
404  (
405  (neg0(iDmdt)*hf2 + pos(iDmdt)*h2)
406  - (pos0(iDmdt)*hf1 + neg(iDmdt)*h1)
407  );
408 
409  volScalarField iDmdtNew(iDmdt);
410 
411  if (phaseChange_)
412  {
413  volScalarField H1(heatTransferModelIter().first()->K(0));
414  volScalarField H2(heatTransferModelIter().second()->K(0));
415 
416  iDmdtNew = (H1*(Tsat - T1) + H2*(Tsat - T2))/L;
417  }
418  else
419  {
420  iDmdtNew == dimensionedScalar(iDmdt.dimensions());
421  }
422 
423  volScalarField H1(heatTransferModelIter().first()->K());
424  volScalarField H2(heatTransferModelIter().second()->K());
425 
426  // Limit the H[12] to avoid /0
427  H1.max(SMALL);
428  H2.max(SMALL);
429 
430  Tf = (H1*T1 + H2*T2 + iDmdtNew*L)/(H1 + H2);
431 
432  Info<< "Tf." << pair.name()
433  << ": min = " << min(Tf.primitiveField())
434  << ", mean = " << average(Tf.primitiveField())
435  << ", max = " << max(Tf.primitiveField())
436  << endl;
437 
438  scalar iDmdtRelax(this->mesh().fieldRelaxationFactor("iDmdt"));
439  iDmdt = (1 - iDmdtRelax)*iDmdt + iDmdtRelax*iDmdtNew;
440 
441  if (phaseChange_)
442  {
443  Info<< "iDmdt." << pair.name()
444  << ": min = " << min(iDmdt.primitiveField())
445  << ", mean = " << average(iDmdt.primitiveField())
446  << ", max = " << max(iDmdt.primitiveField())
447  << ", integral = " << fvc::domainIntegrate(iDmdt).value()
448  << endl;
449  }
450 
451  volScalarField& wDmdt(*this->wDmdt_[pair]);
452  volScalarField& wMDotL(*this->wMDotL_[pair]);
453  wDmdt *= 0.0;
454  wMDotL *= 0.0;
455 
456  bool wallBoilingActive = false;
457 
458  forAllConstIter(phasePair, pair, iter)
459  {
460  const phaseModel& phase = iter();
461  const phaseModel& otherPhase = iter.otherPhase();
462 
463  if
464  (
465  phase.mesh().foundObject<volScalarField>
466  (
467  "alphat." + phase.name()
468  )
469  )
470  {
471  const volScalarField& alphat =
472  phase.mesh().lookupObject<volScalarField>
473  (
474  "alphat." + phase.name()
475  );
476 
477  const fvPatchList& patches = this->mesh().boundary();
478  forAll(patches, patchi)
479  {
480  const fvPatch& currPatch = patches[patchi];
481 
482  if
483  (
484  isA<alphatPhaseChangeWallFunction>
485  (
486  alphat.boundaryField()[patchi]
487  )
488  )
489  {
490  const alphatPhaseChangeWallFunction& PCpatch =
491  refCast<const alphatPhaseChangeWallFunction>
492  (
493  alphat.boundaryField()[patchi]
494  );
495 
496  phasePairKey key(phase.name(), otherPhase.name());
497 
498  if (PCpatch.activePhasePair(key))
499  {
500  wallBoilingActive = true;
501 
502  const scalarField& patchDmdt =
503  PCpatch.dmdt(key);
504  const scalarField& patchMDotL =
505  PCpatch.mDotL(key);
506 
507  const scalar sign
508  (
509  Pair<word>::compare(pair, key)
510  );
511 
512  forAll(patchDmdt, facei)
513  {
514  const label faceCelli =
515  currPatch.faceCells()[facei];
516  wDmdt[faceCelli] -= sign*patchDmdt[facei];
517  wMDotL[faceCelli] -= sign*patchMDotL[facei];
518  }
519  }
520  }
521  }
522  }
523  }
524 
525  if (wallBoilingActive)
526  {
527  Info<< "wDmdt." << pair.name()
528  << ": min = " << min(wDmdt.primitiveField())
529  << ", mean = " << average(wDmdt.primitiveField())
530  << ", max = " << max(wDmdt.primitiveField())
531  << ", integral = " << fvc::domainIntegrate(wDmdt).value()
532  << endl;
533  }
534  }
535 }
536 
537 
538 template<class BasePhaseSystem>
540 {
541  if (BasePhaseSystem::read())
542  {
543  bool readOK = true;
544 
545  // Models ...
546 
547  return readOK;
548  }
549  else
550  {
551  return false;
552  }
553 }
554 
555 
556 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:57
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
L
const vector L(dict.get< vector >("L"))
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:51
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::ThermalPhaseChangePhaseSystem::dmdts
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
Definition: ThermalPhaseChangePhaseSystem.C:194
Foam::ThermalPhaseChangePhaseSystem::iDmdt
tmp< volScalarField > iDmdt(const phasePairKey &key) const
Return the interfacial mass transfer rate for a pair.
Definition: ThermalPhaseChangePhaseSystem.C:38
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::dimDensity
const dimensionSet dimDensity
Foam::saturationModel
Definition: saturationModel.H:53
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:88
Foam::posPart
dimensionedScalar posPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:221
Foam::GeometricField::primitiveField
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Definition: GeometricFieldI.H:53
Foam::neg0
dimensionedScalar neg0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:210
Foam::phaseModel::otherPhase
const phaseModel & otherPhase() const
Return the other phase in this two-phase system.
Foam::compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
Abstract base-class for all alphatWallFunctions supporting phase-change.
Definition: alphatPhaseChangeWallFunctionFvPatchScalarField.H:57
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
he2
volScalarField & he2
Definition: EEqns.H:3
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAllConstIter
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:338
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::ThermalPhaseChangePhaseSystem::read
virtual bool read()
Read base phaseProperties dictionary.
Definition: ThermalPhaseChangePhaseSystem.C:539
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
ThermalPhaseChangePhaseSystem.H
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
alphatPhaseChangeWallFunctionFvPatchScalarField.H
Foam::label
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:62
Foam::Field< scalar >
Foam::fac::average
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition: facAverage.C:47
Foam::ThermalPhaseChangePhaseSystem::~ThermalPhaseChangePhaseSystem
virtual ~ThermalPhaseChangePhaseSystem()
Destructor.
Definition: ThermalPhaseChangePhaseSystem.C:167
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::ThermalPhaseChangePhaseSystem::ThermalPhaseChangePhaseSystem
ThermalPhaseChangePhaseSystem(const fvMesh &)
Construct from fvMesh.
Definition: ThermalPhaseChangePhaseSystem.C:76
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::cellModeller::lookup
const cellModel * lookup(const word &modelName)
Deprecated(2017-11) equivalent to cellModel::ptr static method.
Definition: cellModeller.H:46
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::ThermalPhaseChangePhaseSystem::correctInterfaceThermo
virtual void correctInterfaceThermo()
Correct the interface thermodynamics.
Definition: ThermalPhaseChangePhaseSystem.C:343
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
Foam::phasePairKey
Definition: phasePairKey.H:59
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::blockMeshTools::read
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:33
phase2
phaseModel & phase2
Definition: setRegionFluidFields.H:6
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam::phasePairKey::ordered
bool ordered() const
Return the ordered flag.
Definition: phasePairKey.C:60
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::phase::name
const word & name() const
Definition: phase.H:111
Foam::HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash >
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::phaseModel::name
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:94
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::negPart
dimensionedScalar negPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:232
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash >
Foam::phasePair::name
virtual word name() const
Pair name.
Definition: phasePair.C:92
Foam::ThermalPhaseChangePhaseSystem::saturation
const saturationModel & saturation() const
Return the saturationModel.
Definition: ThermalPhaseChangePhaseSystem.C:175
fvcVolumeIntegrate.H
Volume integrate volField creating a volField.
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::phasePair::phase2
const phaseModel & phase2() const
Return phase 2.
Definition: phasePairI.H:36
Foam::ThermalPhaseChangePhaseSystem::dmdt
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the mass transfer rate for a pair.
Definition: ThermalPhaseChangePhaseSystem.C:184
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:89
Foam::GeometricField::max
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
Definition: GeometricField.C:1112
phase1
phaseModel & phase1
Definition: setRegionFluidFields.H:5
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
massTransfer
phaseSystem::massTransferTable & massTransfer(massTransferPtr())
Foam::phasePair::phase1
const phaseModel & phase1() const
Return phase 1.
Definition: phasePairI.H:30
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
Foam::ThermalPhaseChangePhaseSystem::wDmdt
tmp< volScalarField > wDmdt(const phasePairKey &key) const
Return the boundary mass transfer rate for a pair.
Definition: ThermalPhaseChangePhaseSystem.C:56
Foam::ThermalPhaseChangePhaseSystem::heatTransfer
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
Definition: ThermalPhaseChangePhaseSystem.C:222
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::ThermalPhaseChangePhaseSystem::massTransfer
virtual autoPtr< phaseSystem::massTransferTable > massTransfer() const
Return the mass transfer matrices.
Definition: ThermalPhaseChangePhaseSystem.C:284
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177