momentum.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) 2018-2020 OpenCFD Ltd.
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 
28 #include "momentum.H"
29 #include "fvMesh.H"
30 #include "volFields.H"
31 #include "cellSet.H"
32 #include "cylindricalRotation.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
41  defineTypeNameAndDebug(momentum, 0);
42  addToRunTimeSelectionTable(functionObject, momentum, dictionary);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
48 
49 void Foam::functionObjects::momentum::purgeFields()
50 {
51  obr_.checkOut(scopedName("momentum"));
52  obr_.checkOut(scopedName("angularMomentum"));
53  obr_.checkOut(scopedName("angularVelocity"));
54 }
55 
56 
57 template<class GeoField>
59 Foam::functionObjects::momentum::newField
60 (
61  const word& baseName,
62  const dimensionSet& dims,
63  bool registerObject
64 ) const
65 {
66  return
68  (
69  IOobject
70  (
71  scopedName(baseName),
72  time_.timeName(),
73  mesh_,
76  registerObject
77  ),
78  mesh_,
80  );
81 }
82 
83 
84 void Foam::functionObjects::momentum::calc()
85 {
86  initialise();
87 
88  // Ensure volRegion is properly up-to-date.
89  // Purge old fields if anything has changed (eg, mesh size etc)
90  if (volRegion::update())
91  {
92  purgeFields();
93  }
94 
95  // When field writing is not enabled we need our local storage
96  // for the momentum and angular velocity fields
97  autoPtr<volVectorField> tmomentum, tAngularMom, tAngularVel;
98 
99 
100  // The base fields required
101  const auto& U = lookupObject<volVectorField>(UName_);
102  const auto* rhoPtr = findObject<volScalarField>(rhoName_);
103 
104  const dimensionedScalar rhoRef("rho", dimDensity, rhoRef_);
105 
106  // For quantities such as the mass-averaged angular velocity,
107  // we would calculate the mass per-cell here.
108 
109  // tmp<volScalarField::Internal> tmass =
110  // (
111  // rhoPtr
112  // ? (mesh_.V() * (*rhoPtr))
113  // : (mesh_.V() * rhoRef)
114  // );
115 
116 
117  // Linear momentum
118  // ~~~~~~~~~~~~~~~
119 
120  auto* pmomentum = getObjectPtr<volVectorField>(scopedName("momentum"));
121 
122  if (!pmomentum)
123  {
124  tmomentum = newField<volVectorField>("momentum", dimVelocity*dimMass);
125  pmomentum = tmomentum.get(); // get(), not release()
126  }
127  auto& momentum = *pmomentum;
128 
129  if (rhoPtr)
130  {
131  momentum.ref() = (U * mesh_.V() * (*rhoPtr));
132  }
133  else
134  {
135  momentum.ref() = (U * mesh_.V() * rhoRef);
136  }
137  momentum.correctBoundaryConditions();
138 
139 
140  // Angular momentum
141  // ~~~~~~~~~~~~~~~~
142 
143  auto* pAngularMom =
144  getObjectPtr<volVectorField>(scopedName("angularMomentum"));
145 
146  if (hasCsys_ && !pAngularMom)
147  {
148  tAngularMom =
149  newField<volVectorField>("angularMomentum", dimVelocity*dimMass);
150  pAngularMom = tAngularMom.get(); // get(), not release()
151  }
152  else if (!pAngularMom)
153  {
154  // Angular momentum not requested, but alias to normal momentum
155  // to simplify logic when calculating the summations
156  pAngularMom = pmomentum;
157  }
158  auto& angularMom = *pAngularMom;
159 
160 
161  // Angular velocity
162  // ~~~~~~~~~~~~~~~~
163 
164  auto* pAngularVel =
165  getObjectPtr<volVectorField>(scopedName("angularVelocity"));
166 
167  if (hasCsys_)
168  {
169  if (!pAngularVel)
170  {
171  tAngularVel =
172  newField<volVectorField>("angularVelocity", dimVelocity);
173  pAngularVel = tAngularVel.get(); // get(), not release()
174  }
175  auto& angularVel = *pAngularVel;
176 
177 
178  // Global to local
179 
180  angularVel.primitiveFieldRef() =
181  csys_.invTransform(mesh_.cellCentres(), U.internalField());
182 
183  angularVel.correctBoundaryConditions();
184 
185  if (rhoPtr)
186  {
187  angularMom.ref() = (angularVel * mesh_.V() * (*rhoPtr));
188  }
189  else
190  {
191  angularMom.ref() = (angularVel * mesh_.V() * rhoRef);
192  }
193 
194  angularMom.correctBoundaryConditions();
195  }
196 
197 
198  // Integrate the selection
199 
200  sumMomentum_ = Zero;
201  sumAngularMom_ = Zero;
202 
203  switch (regionType_)
204  {
205  case vrtCellSet:
206  case vrtCellZone:
207  {
208  for (const label celli : cellIDs())
209  {
210  sumMomentum_ += momentum[celli];
211  sumAngularMom_ += angularMom[celli];
212  }
213  break;
214  }
215  case vrtAll:
216  {
217  for (label celli=0; celli < mesh_.nCells(); ++celli)
218  {
219  sumMomentum_ += momentum[celli];
220  sumAngularMom_ += angularMom[celli];
221  }
222  break;
223  }
224  }
225 
226  reduce(sumMomentum_, sumOp<vector>());
227  reduce(sumAngularMom_, sumOp<vector>());
228 }
229 
230 
231 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
232 
234 {
235  if (!writeToFile() || writtenHeader_)
236  {
237  return;
238  }
239 
240  if (hasCsys_)
241  {
242  writeHeader(os, "Momentum, Angular Momentum");
243  writeHeaderValue(os, "origin", csys_.origin());
244  writeHeaderValue(os, "axis", csys_.e3());
245  }
246  else
247  {
248  writeHeader(os, "Momentum");
249  }
250 
251  if (regionType_ != vrtAll)
252  {
254  (
255  os,
256  "Selection " + regionTypeNames_[regionType_]
257  + " = " + regionName_
258  );
259  }
260 
261  writeHeader(os, "");
262  writeCommented(os, "Time");
263  writeTabbed(os, "(momentum_x momentum_y momentum_z)");
264 
265  if (hasCsys_)
266  {
267  writeTabbed(os, "(momentum_r momentum_rtheta momentum_axis)");
268  }
269 
270  writeTabbed(os, "volume");
271  os << endl;
272 
273  writtenHeader_ = true;
274 }
275 
276 
278 {
279  if (initialised_)
280  {
281  return;
282  }
283 
284  if (!foundObject<volVectorField>(UName_))
285  {
287  << "Could not find U: " << UName_ << " in database"
288  << exit(FatalError);
289  }
290 
291 
292  const auto* pPtr = cfindObject<volScalarField>(pName_);
293 
294  if (pPtr && pPtr->dimensions() == dimPressure)
295  {
296  // Compressible - rho is mandatory
297 
298  if (!foundObject<volScalarField>(rhoName_))
299  {
301  << "Could not find rho:" << rhoName_
302  << exit(FatalError);
303  }
304  }
305 
306  initialised_ = true;
307 }
308 
309 
311 {
312  if (log)
313  {
314  Info<< type() << " " << name() << " write:" << nl;
315 
316  Info<< " Sum of Momentum";
317 
318  if (regionType_ != vrtAll)
319  {
320  Info<< ' ' << regionTypeNames_[regionType_]
321  << ' ' << regionName_;
322  }
323 
324  Info<< " (volume " << volRegion::V() << ')' << nl
325  << " linear : " << sumMomentum_ << nl;
326 
327  if (hasCsys_)
328  {
329  Info<< " angular : " << sumAngularMom_ << nl;
330  }
331 
332  Info<< endl;
333  }
334 
335  if (writeToFile())
336  {
337  writeCurrentTime(os);
338 
339  os << tab << sumMomentum_;
340 
341  if (hasCsys_)
342  {
343  os << tab << sumAngularMom_;
344  }
345 
346  os << tab << volRegion::V() << endl;
347  }
348 }
349 
350 
351 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
352 
354 (
355  const word& name,
356  const Time& runTime,
357  const dictionary& dict,
358  bool readFields
359 )
360 :
363  writeFile(mesh_, name, typeName, dict),
364  sumMomentum_(Zero),
365  sumAngularMom_(Zero),
366  UName_(),
367  pName_(),
368  rhoName_(),
369  rhoRef_(1.0),
370  csys_(),
371  hasCsys_(false),
372  writeMomentum_(false),
373  writeVelocity_(false),
374  writePosition_(false),
375  initialised_(false)
376 {
377  if (readFields)
378  {
379  read(dict);
380  Log << endl;
381  }
382 }
383 
384 
386 (
387  const word& name,
388  const objectRegistry& obr,
389  const dictionary& dict,
390  bool readFields
391 )
392 :
395  writeFile(mesh_, name, typeName, dict),
396  sumMomentum_(Zero),
397  sumAngularMom_(Zero),
398  UName_(),
399  pName_(),
400  rhoName_(),
401  rhoRef_(1.0),
402  csys_(),
403  hasCsys_(false),
404  writeMomentum_(false),
405  writeVelocity_(false),
406  writePosition_(false),
407  initialised_(false)
408 {
409  if (readFields)
410  {
411  read(dict);
412  Log << endl;
413  }
414 }
415 
416 
417 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
418 
420 {
424 
425  initialised_ = false;
426 
427  Info<< type() << " " << name() << ":" << nl;
428 
429  // Optional entries U and p
430  UName_ = dict.getOrDefault<word>("U", "U");
431  pName_ = dict.getOrDefault<word>("p", "p");
432  rhoName_ = dict.getOrDefault<word>("rho", "rho");
433  rhoRef_ = dict.getOrDefault<scalar>("rhoRef", 1.0);
434  hasCsys_ = dict.getOrDefault("cylindrical", false);
435 
436  if (hasCsys_)
437  {
439  }
440 
441  writeMomentum_ = dict.getOrDefault("writeMomentum", false);
442  writeVelocity_ = dict.getOrDefault("writeVelocity", false);
443  writePosition_ = dict.getOrDefault("writePosition", false);
444 
445  Info<<"Integrating for selection: "
446  << regionTypeNames_[regionType_]
447  << " (" << regionName_ << ")" << nl;
448 
449  if (writeMomentum_)
450  {
451  Info<< " Momentum fields will be written" << endl;
452 
453  mesh_.objectRegistry::store
454  (
455  newField<volVectorField>("momentum", dimVelocity*dimMass)
456  );
457 
458  if (hasCsys_)
459  {
460  mesh_.objectRegistry::store
461  (
462  newField<volVectorField>("angularMomentum", dimVelocity*dimMass)
463  );
464  }
465  }
466 
467  if (hasCsys_)
468  {
469  if (writeVelocity_)
470  {
471  Info<< " Angular velocity will be written" << endl;
472 
473  mesh_.objectRegistry::store
474  (
475  newField<volVectorField>("angularVelocity", dimVelocity)
476  );
477  }
478 
479  if (writePosition_)
480  {
481  Info<< " Angular position will be written" << endl;
482  }
483  }
484 
485  return true;
486 }
487 
488 
490 {
491  calc();
492 
493  if (Pstream::master())
494  {
495  writeFileHeader(file());
496 
497  writeValues(file());
498 
499  Log << endl;
500  }
501 
502  // Write state/results information
503  setResult("momentum_x", sumMomentum_[0]);
504  setResult("momentum_y", sumMomentum_[1]);
505  setResult("momentum_z", sumMomentum_[2]);
506 
507  setResult("momentum_r", sumAngularMom_[0]);
508  setResult("momentum_rtheta", sumAngularMom_[1]);
509  setResult("momentum_axis", sumAngularMom_[2]);
510 
511  return true;
512 }
513 
514 
516 {
517  if (writeMomentum_ || (hasCsys_ && (writeVelocity_ || writePosition_)))
518  {
519  Log << "Writing fields" << nl;
520 
521  const volVectorField* fieldPtr;
522 
523  fieldPtr = findObject<volVectorField>(scopedName("momentum"));
524  if (fieldPtr) fieldPtr->write();
525 
526  fieldPtr = findObject<volVectorField>(scopedName("angularMomentum"));
527  if (fieldPtr) fieldPtr->write();
528 
529  fieldPtr = findObject<volVectorField>(scopedName("angularVelocity"));
530  if (fieldPtr) fieldPtr->write();
531 
532  if (hasCsys_ && writePosition_)
533  {
534  // Clunky, but currently no simple means of handling
535  // component-wise conversion and output
536 
537  auto cyl_r = newField<volScalarField>("cyl_r", dimLength);
538  auto cyl_t = newField<volScalarField>("cyl_theta", dimless);
539  auto cyl_z = newField<volScalarField>("cyl_z", dimLength);
540 
541  // Internal
542  {
543  const auto& pts = mesh_.cellCentres();
544  const label len = pts.size();
545 
546  UList<scalar>& r = cyl_r->primitiveFieldRef(false);
547  UList<scalar>& t = cyl_t->primitiveFieldRef(false);
548  UList<scalar>& z = cyl_z->primitiveFieldRef(false);
549 
550  for (label i=0; i < len; ++i)
551  {
552  point p(csys_.localPosition(pts[i]));
553 
554  r[i] = p.x();
555  t[i] = p.y();
556  z[i] = p.z();
557  }
558  }
559 
560  // Boundary
561  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
562 
563  forAll(pbm, patchi)
564  {
565  const auto& pts = pbm[patchi].faceCentres();
566  const label len = pts.size();
567 
568  UList<scalar>& r = cyl_r->boundaryFieldRef(false)[patchi];
569  UList<scalar>& t = cyl_t->boundaryFieldRef(false)[patchi];
570  UList<scalar>& z = cyl_z->boundaryFieldRef(false)[patchi];
571 
572  for (label i=0; i < len; ++i)
573  {
574  point p(csys_.localPosition(pts[i]));
575 
576  r[i] = p.x();
577  t[i] = p.y();
578  z[i] = p.z();
579  }
580  }
581 
582  cyl_r->write();
583  cyl_t->write();
584  cyl_z->write();
585  }
586  }
587 
588  return true;
589 }
590 
591 
593 {
595  purgeFields(); // Mesh motion makes calculated fields dubious
596 }
597 
598 
600 {
602  purgeFields(); // Mesh motion makes calculated fields dubious
603 }
604 
605 
606 // ************************************************************************* //
Foam::functionObjects::momentum::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: momentum.C:592
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
volFields.H
Foam::functionObjects::momentum::read
virtual bool read(const dictionary &)
Read the momentum data.
Definition: momentum.C:419
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::dimPressure
const dimensionSet dimPressure
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::functionObjects::regionFunctionObject::obr_
const objectRegistry & obr_
Reference to the region objectRegistry.
Definition: regionFunctionObject.H:103
Foam::functionObjects::momentum::write
virtual bool write()
Write the momentum, possibly angular momentum and velocity.
Definition: momentum.C:515
p
volScalarField & p
Definition: createFieldRefs.H:8
Log
#define Log
Definition: PDRblock.C:35
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::functionObjects::volRegion::update
bool update()
Update the cached values as required.
Definition: volRegion.C:228
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::functionObjects::volRegion
Volume (cell) region selection class.
Definition: volRegion.H:115
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::functionObjects::momentum::scopedName
word scopedName(const word &base) const
Return a scoped name.
Definition: momentum.H:331
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:458
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::functionObjects::momentum::momentum
momentum(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from Time and dictionary.
Definition: momentum.C:354
Foam::writeHeader
static void writeHeader(Ostream &os, const word &fieldName)
Definition: rawSurfaceWriterImpl.C:49
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::functionObjects::volRegion::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: volRegion.C:240
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
rhoPtr
Info<< "Reading mechanical properties\n"<< endl;IOdictionary mechanicalProperties(IOobject("mechanicalProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));const dictionary &rhoDict(mechanicalProperties.subDict("rho"));word rhoType(rhoDict.get< word >"type"));autoPtr< volScalarField > rhoPtr
Definition: readMechanicalProperties.H:18
Foam::functionObjects::volRegion::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: volRegion.C:167
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:214
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::functionObjects::momentum::execute
virtual bool execute()
Calculate and report the integral momentum.
Definition: momentum.C:489
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::autoPtr::get
T * get() noexcept
Return pointer to managed object without nullptr checking.
Definition: autoPtr.H:158
Foam::coordSystem::cylindrical
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
Definition: cylindricalCS.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::functionObjects::momentum::writeFileHeader
virtual void writeFileHeader(Ostream &os)
Output file header information.
Definition: momentum.C:233
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::functionObjects::momentum::movePoints
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: momentum.C:599
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
cylindricalRotation.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::tab
constexpr char tab
Definition: Ostream.H:384
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::functionObjects::momentum::initialise
void initialise()
Initialise the fields.
Definition: momentum.C:277
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::Vector< scalar >
Foam::functionObjects::volRegion::movePoints
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: volRegion.C:246
Foam::functionObjects::momentum::writeValues
void writeValues(Ostream &os)
Write momentum data.
Definition: momentum.C:310
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Reference to the fvMesh.
Definition: fvMeshFunctionObject.H:73
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::functionObjects::readFields
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:155
Foam::objectRegistry::checkOut
bool checkOut(regIOobject *io) const
Definition: objectRegistry.C:271
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::functionObjects::log
Computes the natural logarithm of an input volScalarField.
Definition: log.H:227
Foam::functionObjects::writeFile
Base class for writing single files from the function objects.
Definition: writeFile.H:119
momentum.H
cellSet.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::functionObjects::volRegion::V
scalar V() const
Return total volume of the selected region.
Definition: volRegionI.H:62