sampledCuttingPlane.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "sampledCuttingPlane.H"
30 #include "dictionary.H"
31 #include "volFields.H"
32 #include "volPointInterpolation.H"
34 #include "fvMesh.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(sampledCuttingPlane, 0);
42  (
43  sampledSurface,
44  sampledCuttingPlane,
45  word,
46  cuttingPlane
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::sampledCuttingPlane::checkBoundsIntersection
54 (
55  const plane& pln,
56  const boundBox& meshBb
57 ) const
58 {
59  // Verify specified bounding box
60  if (bounds_.valid())
61  {
62  // Bounding box does not overlap with (global) mesh!
63  if (!bounds_.overlaps(meshBb))
64  {
66  << nl
67  << name() << " : "
68  << "Bounds " << bounds_
69  << " do not overlap the mesh bounding box " << meshBb
70  << nl << endl;
71  }
72 
73  // Plane does not intersect the bounding box
74  if (!bounds_.intersects(pln))
75  {
77  << nl
78  << name() << " : "
79  << "Plane "<< pln << " does not intersect the bounds "
80  << bounds_
81  << nl << endl;
82  }
83  }
84 
85  // Plane does not intersect the (global) mesh!
86  if (!meshBb.intersects(pln))
87  {
89  << nl
90  << name() << " : "
91  << "Plane "<< pln << " does not intersect the mesh bounds "
92  << meshBb
93  << nl << endl;
94  }
95 }
96 
97 
98 void Foam::sampledCuttingPlane::createGeometry()
99 {
100  if (debug)
101  {
102  Pout<< "sampledCuttingPlane::createGeometry :updating geometry."
103  << endl;
104  }
105 
106  // Clear any stored topologies
107  isoSurfPtr_.clear();
108  isoSurfCellPtr_.clear();
109  isoSurfTopoPtr_.clear();
110  pointDistance_.clear();
111  cellDistancePtr_.clear();
112 
113  // Clear derived data
114  clearGeom();
115 
116  const fvMesh& fvm = static_cast<const fvMesh&>(this->mesh());
117 
118  // Get sub-mesh if any
119  if
120  (
121  (-1 != mesh().cellZones().findIndex(zoneNames_))
122  && subMeshPtr_.empty()
123  )
124  {
125  const polyBoundaryMesh& patches = mesh().boundaryMesh();
126 
127  // Patch to put exposed internal faces into
128  const label exposedPatchi = patches.findPatchID(exposedPatchName_);
129 
130  bitSet cellsToSelect = mesh().cellZones().selection(zoneNames_);
131 
132  DebugInfo
133  << "Allocating subset of size "
134  << cellsToSelect.count()
135  << " with exposed faces into patch "
136  << patches[exposedPatchi].name() << endl;
137 
138 
139  // If we will use a fvMeshSubset so can apply bounds as well to make
140  // the initial selection smaller.
141  if (bounds_.valid() && cellsToSelect.any())
142  {
143  const auto& cellCentres = fvm.C();
144 
145  for (const label celli : cellsToSelect)
146  {
147  const point& cc = cellCentres[celli];
148 
149  if (!bounds_.contains(cc))
150  {
151  cellsToSelect.unset(celli);
152  }
153  }
154 
155  DebugInfo
156  << "Bounded subset of size "
157  << cellsToSelect.count() << endl;
158  }
159 
160  subMeshPtr_.reset(new fvMeshSubset(fvm, cellsToSelect, exposedPatchi));
161  }
162 
163 
164  // Select either the submesh or the underlying mesh
165  const fvMesh& mesh =
166  (
167  subMeshPtr_.valid()
168  ? subMeshPtr_().subMesh()
169  : fvm
170  );
171 
172  checkBoundsIntersection(plane_, mesh.bounds());
173 
174 
175  // Distance to cell centres
176  // ~~~~~~~~~~~~~~~~~~~~~~~~
177 
178  cellDistancePtr_.reset
179  (
180  new volScalarField
181  (
182  IOobject
183  (
184  "cellDistance",
185  mesh.time().timeName(),
186  mesh.time(),
189  false
190  ),
191  mesh,
193  )
194  );
195  volScalarField& cellDistance = cellDistancePtr_();
196 
197  // Internal field
198  {
199  const auto& cc = mesh.cellCentres();
200  scalarField& fld = cellDistance.primitiveFieldRef();
201 
202  forAll(cc, i)
203  {
204  fld[i] = plane_.signedDistance(cc[i]);
205  }
206  }
207 
208  volScalarField::Boundary& cellDistanceBf =
209  cellDistance.boundaryFieldRef();
210 
211  // Patch fields
212  {
213  forAll(cellDistanceBf, patchi)
214  {
215  if
216  (
217  isA<emptyFvPatchScalarField>
218  (
219  cellDistanceBf[patchi]
220  )
221  )
222  {
223  cellDistanceBf.set
224  (
225  patchi,
226  new calculatedFvPatchScalarField
227  (
228  mesh.boundary()[patchi],
229  cellDistance
230  )
231  );
232 
233  const polyPatch& pp = mesh.boundary()[patchi].patch();
234  pointField::subField cc = pp.patchSlice(mesh.faceCentres());
235 
236  fvPatchScalarField& fld = cellDistanceBf[patchi];
237  fld.setSize(pp.size());
238  forAll(fld, i)
239  {
240  fld[i] = plane_.signedDistance(cc[i]);
241  }
242  }
243  else
244  {
245  // Other side cell centres?
246  const pointField& cc = mesh.C().boundaryField()[patchi];
247  fvPatchScalarField& fld = cellDistanceBf[patchi];
248 
249  forAll(fld, i)
250  {
251  fld[i] = plane_.signedDistance(cc[i]);
252  }
253  }
254  }
255  }
256 
257 
258  // On processor patches the mesh.C() will already be the cell centre
259  // on the opposite side so no need to swap cellDistance.
260 
261 
262  // Distance to points
263  pointDistance_.setSize(mesh.nPoints());
264  {
265  const pointField& pts = mesh.points();
266 
267  forAll(pointDistance_, i)
268  {
269  pointDistance_[i] = plane_.signedDistance(pts[i]);
270  }
271  }
272 
273 
274  if (debug)
275  {
276  Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
277  cellDistance.write();
278  pointScalarField pDist
279  (
280  IOobject
281  (
282  "pointDistance",
283  mesh.time().timeName(),
284  mesh.time(),
287  false
288  ),
291  );
292  pDist.primitiveFieldRef() = pointDistance_;
293 
294  Pout<< "Writing point distance:" << pDist.objectPath() << endl;
295  pDist.write();
296  }
297 
298 
299  // Direct from cell field and point field.
300  if (isoAlgo_ == isoSurfaceBase::ALGO_CELL)
301  {
302  isoSurfCellPtr_.reset
303  (
304  new isoSurfaceCell
305  (
306  fvm,
307  cellDistance,
308  pointDistance_,
309  0,
310  filter_,
311  bounds_,
312  mergeTol_
313  )
314  );
315  }
316  else if (isoAlgo_ == isoSurfaceBase::ALGO_TOPO)
317  {
318  isoSurfTopoPtr_.reset
319  (
320  new isoSurfaceTopo
321  (
322  fvm,
323  cellDistance,
324  pointDistance_,
325  0,
326  filter_,
327  bounds_
328  )
329  );
330  }
331  else
332  {
333  isoSurfPtr_.reset
334  (
335  new isoSurface
336  (
337  cellDistance,
338  pointDistance_,
339  0,
340  filter_,
341  bounds_,
342  mergeTol_
343  )
344  );
345  }
346 
347  if (debug)
348  {
349  print(Pout);
350  Pout<< endl;
351  }
352 }
353 
354 
355 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
356 
358 (
359  const word& name,
360  const polyMesh& mesh,
361  const dictionary& dict
362 )
363 :
365  plane_(dict),
366  bounds_(dict.getOrDefault("bounds", boundBox::invertedBox)),
367  mergeTol_(dict.getOrDefault<scalar>("mergeTol", 1e-6)),
368  isoAlgo_
369  (
370  isoSurfaceBase::algorithmNames.getOrDefault
371  (
372  "isoAlgorithm",
373  dict,
375  )
376  ),
377  filter_
378  (
380  (
381  dict,
383  )
384  ),
385  average_(dict.getOrDefault("average", false)),
386  zoneNames_(),
387  exposedPatchName_(),
388  needsUpdate_(true),
389  subMeshPtr_(nullptr),
390  cellDistancePtr_(nullptr),
391  isoSurfPtr_(nullptr),
392  isoSurfCellPtr_(nullptr),
393  isoSurfTopoPtr_(nullptr)
394 {
395  if (!dict.readIfPresent("zones", zoneNames_) && dict.found("zone"))
396  {
397  zoneNames_.resize(1);
398  dict.readEntry("zone", zoneNames_.first());
399  }
400 
401  if (-1 != mesh.cellZones().findIndex(zoneNames_))
402  {
403  dict.readEntry("exposedPatchName", exposedPatchName_);
404 
405  if (-1 == mesh.boundaryMesh().findPatchID(exposedPatchName_))
406  {
408  << "Cannot find patch " << exposedPatchName_
409  << " in which to put exposed faces." << endl
410  << "Valid patches are " << mesh.boundaryMesh().names()
411  << exit(FatalError);
412  }
413 
414  DebugInfo
415  << "Restricting to cellZone(s) " << flatOutput(zoneNames_)
416  << " with exposed internal faces into patch "
417  << exposedPatchName_ << endl;
418  }
419 }
420 
421 
422 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
423 
425 {
426  return needsUpdate_;
427 }
428 
429 
431 {
432  if (debug)
433  {
434  Pout<< "sampledCuttingPlane::expire :"
435  << " needsUpdate:" << needsUpdate_ << endl;
436  }
437 
438  // Clear derived data
439  clearGeom();
440 
441  // Already marked as expired
442  if (needsUpdate_)
443  {
444  return false;
445  }
446 
447  needsUpdate_ = true;
448  return true;
449 }
450 
451 
453 {
454  if (debug)
455  {
456  Pout<< "sampledCuttingPlane::update :"
457  << " needsUpdate:" << needsUpdate_ << endl;
458  }
459 
460  if (!needsUpdate_)
461  {
462  return false;
463  }
464 
465  createGeometry();
466 
467  needsUpdate_ = false;
468  return true;
469 }
470 
471 
474 (
475  const interpolation<scalar>& sampler
476 ) const
477 {
478  return sampleOnFaces(sampler);
479 }
480 
481 
484 (
485  const interpolation<vector>& sampler
486 ) const
487 {
488  return sampleOnFaces(sampler);
489 }
490 
491 
494 (
495  const interpolation<sphericalTensor>& sampler
496 ) const
497 {
498  return sampleOnFaces(sampler);
499 }
500 
501 
504 (
505  const interpolation<symmTensor>& sampler
506 ) const
507 {
508  return sampleOnFaces(sampler);
509 }
510 
511 
514 (
515  const interpolation<tensor>& sampler
516 ) const
517 {
518  return sampleOnFaces(sampler);
519 }
520 
521 
524 (
525  const interpolation<scalar>& interpolator
526 ) const
527 {
528  return sampleOnPoints(interpolator);
529 }
530 
531 
534 (
535  const interpolation<vector>& interpolator
536 ) const
537 {
538  return sampleOnPoints(interpolator);
539 }
540 
541 
544 (
545  const interpolation<sphericalTensor>& interpolator
546 ) const
547 {
548  return sampleOnPoints(interpolator);
549 }
550 
551 
554 (
555  const interpolation<symmTensor>& interpolator
556 ) const
557 {
558  return sampleOnPoints(interpolator);
559 }
560 
561 
564 (
565  const interpolation<tensor>& interpolator
566 ) const
567 {
568  return sampleOnPoints(interpolator);
569 }
570 
571 
573 {
574  os << "sampledCuttingPlane: " << name() << " :"
575  << " plane:" << plane_
576  << " faces:" << faces().size()
577  << " points:" << points().size();
578 }
579 
580 
581 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:40
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
volFields.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:46
Foam::isoSurfaceBase::filterType::DIAGCELL
Remove points from pyramid edges and face-diagonals.
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:359
Foam::boundBox::invertedBox
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Foam::isoSurfaceBase::ALGO_CELL
Definition: isoSurfaceBase.H:81
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::isoSurfaceBase::ALGO_POINT
Definition: isoSurfaceBase.H:80
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:483
Foam::isoSurfaceBase::algorithmNames
static const Enum< algorithmType > algorithmNames
Names for the iso-surface algorithms.
Definition: isoSurfaceBase.H:100
Foam::ZoneMesh::findIndex
label findIndex(const keyType &key) const
Return zone index for the first match, return -1 if not found.
Definition: ZoneMesh.C:449
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
meshBb
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:593
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::sampledCuttingPlane::print
virtual void print(Ostream &os) const
Print information.
Definition: sampledCuttingPlane.C:572
Foam::boundBox::overlaps
bool overlaps(const boundBox &bb) const
Overlaps/touches boundingBox?
Definition: boundBoxI.H:221
sampledCuttingPlane.H
Foam::sampledCuttingPlane::needsUpdate
virtual bool needsUpdate() const
Does the surface need an update?
Definition: sampledCuttingPlane.C:424
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::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::addNamedToRunTimeSelectionTable
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::sampledCuttingPlane::update
virtual bool update()
Update the surface as required.
Definition: sampledCuttingPlane.C:452
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:120
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:766
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:358
Foam::interpolation< scalar >
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
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::pointScalarField
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:51
Foam::sampledCuttingPlane::sampledCuttingPlane
sampledCuttingPlane(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Definition: sampledCuttingPlane.C:358
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:735
volPointInterpolation.H
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:441
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:176
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::flatOutput
FlatOutput< Container > flatOutput(const Container &obj, label len=0)
Global flatOutput function.
Definition: FlatOutput.H:85
Foam::boundBox::intersects
bool intersects(const plane &pln) const
Does plane intersect this bounding box.
Definition: boundBox.C:200
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:145
Foam::sampledSurface::name
const word & name() const
Name of surface.
Definition: sampledSurface.H:308
Foam::BitOps::print
Ostream & print(Ostream &os, UIntType value, char off='0', char on='1')
Print 0/1 bits in the (unsigned) integral type.
Definition: BitOps.H:172
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::primitiveMesh::nPoints
label nPoints() const
Number of mesh points.
Definition: primitiveMeshI.H:37
dictionary.H
Foam::isoSurfaceBase::getFilterType
static filterType getFilterType(const dictionary &dict, const isoSurfaceBase::filterType deflt)
Get 'regularise' as bool or enumeration.
Definition: isoSurfaceBase.C:62
Foam::primitiveMesh::reset
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
Definition: primitiveMesh.C:208
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::Field::subField
SubField< Type > subField
Declare type of subField.
Definition: Field.H:90
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::sampledCuttingPlane::sample
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
Sample volume field onto surface faces.
Definition: sampledCuttingPlane.C:474
Foam::ZoneMesh::selection
bitSet selection(const labelUList &zoneIds) const
Definition: ZoneMesh.C:528
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::isoSurfaceBase::ALGO_TOPO
Definition: isoSurfaceBase.H:82
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::boundBox::valid
bool valid() const
Bounding box is non-inverted.
Definition: boundBoxI.H:76
Foam::sampledCuttingPlane::expire
virtual bool expire()
Mark the surface as needing an update.
Definition: sampledCuttingPlane.C:430
Foam::sampledSurface::interpolate
bool interpolate() const
Interpolation to nodes requested for surface.
Definition: sampledSurface.H:326
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62