setFields.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) 2022 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 Application
28  setFields
29 
30 Group
31  grpPreProcessingUtilities
32 
33 Description
34  Set values on a selected set of cells/patch-faces via a dictionary.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "argList.H"
39 #include "Time.H"
40 #include "fvMesh.H"
41 #include "faMesh.H"
42 #include "topoSetSource.H"
43 #include "cellSet.H"
44 #include "faceSet.H"
45 #include "volFields.H"
46 #include "areaFields.H"
47 
48 using namespace Foam;
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 //- Simple tuple of field type and name (read from stream)
53 class fieldDescription
54 {
55  word type_;
56  word name_;
57 
58 public:
59 
60  const word& type() const noexcept { return type_; }
61  const word& name() const noexcept { return name_; }
62 
63  explicit fieldDescription(Istream& is)
64  {
65  is >> type_;
66  is >> name_;
67 
68  // Eg, read as "volScalarFieldValue", but change to "volScalarField"
69  if (type_.ends_with("Value"))
70  {
71  type_.erase(type_.size()-5);
72  }
73  }
74 };
75 
76 
77 // Consume unused field information
78 template<class Type>
79 bool consumeUnusedType(const fieldDescription& fieldDesc, Istream& is)
80 {
83  //? typedef GeometricField<Type, faePatchField, areaMesh> fieldType3;
84  //? typedef GeometricField<Type, fvsPatchField, volMesh> fieldType4;
85 
86  if
87  (
88  fieldDesc.type() == fieldType1::typeName
89  || fieldDesc.type() == fieldType2::typeName
90  )
91  {
92  (void) pTraits<Type>(is);
93  return true;
94  }
95 
96  return false;
97 }
98 
99 
100 // Consume unused field information
101 static bool consumeUnused(const fieldDescription& fieldDesc, Istream& is)
102 {
103  return
104  (
105  consumeUnusedType<scalar>(fieldDesc, is)
106  || consumeUnusedType<vector>(fieldDesc, is)
107  || consumeUnusedType<sphericalTensor>(fieldDesc, is)
108  || consumeUnusedType<symmTensor>(fieldDesc, is)
109  || consumeUnusedType<tensor>(fieldDesc, is)
110  );
111 }
112 
113 
114 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
115 
116 // Setting volume fields
117 template<class Type>
118 bool setCellFieldType
119 (
120  const fieldDescription& fieldDesc,
121  const fvMesh& mesh,
122  const labelList& selectedCells,
123  Istream& is
124 )
125 {
127 
128  if (fieldDesc.type() != fieldType::typeName)
129  {
130  return false;
131  }
132 
133  // Get value from stream
134  const Type fieldValue = pTraits<Type>(is);
135 
136 
137  // Check the current time directory
138  IOobject fieldHeader
139  (
140  fieldDesc.name(),
141  mesh.thisDb().time().timeName(),
142  mesh.thisDb(),
144  );
145 
146  bool found = fieldHeader.typeHeaderOk<fieldType>(true);
147 
148  if (!found)
149  {
150  // Fallback to "constant" directory
151  fieldHeader = IOobject
152  (
153  fieldDesc.name(),
154  mesh.thisDb().time().constant(),
155  mesh.thisDb(),
157  );
158  found = fieldHeader.typeHeaderOk<fieldType>(true);
159  }
160 
161  // Field exists
162  if (found)
163  {
164  Info<< " - set internal values of "
165  << fieldHeader.headerClassName()
166  << ": " << fieldDesc.name()
167  << " = " << fieldValue << endl;
168 
169  fieldType field(fieldHeader, mesh, false);
170 
171  if (isNull(selectedCells) || selectedCells.size() == field.size())
172  {
173  field.primitiveFieldRef() = fieldValue;
174  }
175  else
176  {
177  for (const label celli : selectedCells)
178  {
179  field[celli] = fieldValue;
180  }
181  }
182 
183  // Make boundary fields consistent - treat like zeroGradient
184  for (auto& pfld : field.boundaryFieldRef())
185  {
186  pfld = pfld.patchInternalField();
187  }
188 
189  if (!field.write())
190  {
192  << "Failed writing field " << field.name() << endl;
193  }
194  }
195  else
196  {
197  Warning
198  << "Field " << fieldDesc.name() << " not found" << endl;
199  }
200 
201  return true;
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 // Setting finite-area face fields
208 template<class Type>
209 bool setAreaFieldType
210 (
211  const fieldDescription& fieldDesc,
212  const faMesh& mesh,
213  const labelList& selectedFaces,
214  Istream& is
215 )
216 {
218 
219  if (fieldDesc.type() != fieldType::typeName)
220  {
221  return false;
222  }
223 
224  // Get value from stream
225  const Type fieldValue = pTraits<Type>(is);
226 
227 
228  // Check the current time directory
229  IOobject fieldHeader
230  (
231  fieldDesc.name(),
232  mesh.thisDb().time().timeName(),
233  mesh.thisDb(),
235  );
236 
237  bool found = fieldHeader.typeHeaderOk<fieldType>(true);
238 
239  if (!found)
240  {
241  // Fallback to "constant" directory
242  fieldHeader = IOobject
243  (
244  fieldDesc.name(),
245  mesh.thisDb().time().constant(),
246  mesh.thisDb(),
248  );
249  found = fieldHeader.typeHeaderOk<fieldType>(true);
250  }
251 
252  // Field exists
253  if (found)
254  {
255  Info<< " - set internal values of "
256  << fieldHeader.headerClassName()
257  << ": " << fieldDesc.name()
258  << " = " << fieldValue << endl;
259 
260  fieldType field(fieldHeader, mesh);
261 
262  if (isNull(selectedFaces) || selectedFaces.size() == field.size())
263  {
264  field.primitiveFieldRef() = fieldValue;
265  }
266  else
267  {
268  for (const label facei : selectedFaces)
269  {
270  field[facei] = fieldValue;
271  }
272  }
273 
274  if (!field.write())
275  {
277  << "Failed writing field " << field.name() << endl;
278  }
279  }
280  else
281  {
282  Warning
283  << "Field " << fieldDesc.name() << " not found" << endl;
284  }
285 
286  return true;
287 }
288 
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 // Setting volume boundary fields
293 template<class Type>
294 bool setFaceFieldType
295 (
296  const fieldDescription& fieldDesc,
297  const fvMesh& mesh,
298  const labelList& selectedFaces,
299  Istream& is
300 )
301 {
303 
304  if (fieldDesc.type() != fieldType::typeName)
305  {
306  return false;
307  }
308 
309  // Get value from stream
310  const Type fieldValue = pTraits<Type>(is);
311 
312 
313  // Check the current time directory
314  IOobject fieldHeader
315  (
316  fieldDesc.name(),
317  mesh.thisDb().time().timeName(),
318  mesh.thisDb(),
320  );
321 
322  bool found = fieldHeader.typeHeaderOk<fieldType>(true);
323 
324  if (!found)
325  {
326  // Fallback to "constant" directory
327  fieldHeader = IOobject
328  (
329  fieldDesc.name(),
330  mesh.thisDb().time().constant(),
331  mesh.thisDb(),
333  );
334  found = fieldHeader.typeHeaderOk<fieldType>(true);
335  }
336 
337  // Field exists
338  if (found)
339  {
340  Info<< " - set boundary values of "
341  << fieldHeader.headerClassName()
342  << ": " << fieldDesc.name()
343  << " = " << fieldValue << endl;
344 
345  fieldType field(fieldHeader, mesh);
346 
347  // Create flat list of selected faces and their value.
348  Field<Type> allBoundaryValues(mesh.nBoundaryFaces());
349  forAll(field.boundaryField(), patchi)
350  {
352  (
353  allBoundaryValues,
354  field.boundaryField()[patchi].size(),
355  field.boundaryField()[patchi].patch().start()
356  - mesh.nInternalFaces()
357  ) = field.boundaryField()[patchi];
358  }
359 
360  // Override
361  unsigned hasWarned = 0;
362  labelList nChanged
363  (
364  returnReduce(field.boundaryField().size(), maxOp<label>()),
365  Zero
366  );
367 
368  for (const label facei : selectedFaces)
369  {
370  const label bFacei = facei-mesh.nInternalFaces();
371 
372  if (bFacei < 0)
373  {
374  if (!(hasWarned & 1))
375  {
376  hasWarned |= 1;
378  << "Ignoring internal face " << facei
379  << ". Suppressing further warnings." << endl;
380  }
381  }
382  else if (bFacei >= mesh.nBoundaryFaces())
383  {
384  if (!(hasWarned & 2))
385  {
386  hasWarned |= 2;
388  << "Ignoring out-of-range face " << facei
389  << ". Suppressing further warnings." << endl;
390  }
391  }
392  else
393  {
394  label patchi = mesh.boundaryMesh().patchID()[bFacei];
395 
396  allBoundaryValues[bFacei] = fieldValue;
397  ++nChanged[patchi];
398  }
399  }
400 
402 
403  auto& fieldBf = field.boundaryFieldRef();
404 
405  // Reassign.
406  forAll(field.boundaryField(), patchi)
407  {
408  if (nChanged[patchi] > 0)
409  {
410  Info<< " On patch "
411  << field.boundaryField()[patchi].patch().name()
412  << " set " << nChanged[patchi] << " values" << endl;
413 
414  fieldBf[patchi] == SubField<Type>
415  (
416  allBoundaryValues,
417  fieldBf[patchi].size(),
418  fieldBf[patchi].patch().start()
419  - mesh.nInternalFaces()
420  );
421  }
422  }
423 
424  if (!field.write())
425  {
427  << "Failed writing field " << field.name() << endl;
428  }
429  }
430  else
431  {
432  Warning
433  << "Field " << fieldDesc.name() << " not found" << endl;
434  }
435 
436  return true;
437 }
438 
439 
440 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
441 
442 // Dispatcher for setting volume fields
443 struct setCellField
444 {
445  autoPtr<setCellField> clone() const { return nullptr; } // placeholder
446 
447  static bool apply
448  (
449  const fieldDescription& fieldDesc,
450  const fvMesh& m,
451  const labelList& selectedCells,
452  Istream& is
453  )
454  {
455  return
456  (
457  setCellFieldType<scalar>(fieldDesc, m, selectedCells, is)
458  || setCellFieldType<vector>(fieldDesc, m, selectedCells, is)
459  || setCellFieldType<sphericalTensor>(fieldDesc, m, selectedCells, is)
460  || setCellFieldType<symmTensor>(fieldDesc, m, selectedCells, is)
461  || setCellFieldType<tensor>(fieldDesc, m, selectedCells, is)
462  );
463  }
464 
465  class iNew
466  {
467  const fvMesh& mesh_;
468  const labelList& selected_;
469 
470  public:
471 
472  iNew(const fvMesh& mesh, const labelList& selectedCells)
473  :
474  mesh_(mesh),
475  selected_(selectedCells)
476  {}
477 
478  autoPtr<setCellField> operator()(Istream& is) const
479  {
480  const fieldDescription fieldDesc(is);
481 
482  bool ok = setCellField::apply(fieldDesc, mesh_, selected_, is);
483 
484  if (!ok)
485  {
486  ok = consumeUnused(fieldDesc, is);
487 
488  if (ok)
489  {
490  // Not meant for us
491  Info<< "Skip " << fieldDesc.type()
492  << " for finite-volume" << nl;
493  }
494  else
495  {
497  << "Unsupported field type: "
498  << fieldDesc.type() << endl;
499  }
500  }
501 
502  return nullptr; // Irrelevant return value
503  }
504  };
505 };
506 
507 
508 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
509 
510 // Dispatcher for setting volume boundary fields
511 struct setFaceField
512 {
513  autoPtr<setFaceField> clone() const { return nullptr; } // placeholder
514 
515  static bool apply
516  (
517  const fieldDescription& fieldDesc,
518  const fvMesh& m,
519  const labelList& selectedFaces,
520  Istream& is
521  )
522  {
523  return
524  (
525  setFaceFieldType<scalar>(fieldDesc, m, selectedFaces, is)
526  || setFaceFieldType<vector>(fieldDesc, m, selectedFaces, is)
527  || setFaceFieldType<sphericalTensor>(fieldDesc, m, selectedFaces, is)
528  || setFaceFieldType<symmTensor>(fieldDesc, m, selectedFaces, is)
529  || setFaceFieldType<tensor>(fieldDesc, m, selectedFaces, is)
530  );
531  }
532 
533  class iNew
534  {
535  const fvMesh& mesh_;
536  const labelList& selected_;
537 
538  public:
539 
540  iNew(const fvMesh& mesh, const labelList& selectedFaces)
541  :
542  mesh_(mesh),
543  selected_(selectedFaces)
544  {}
545 
546  autoPtr<setFaceField> operator()(Istream& is) const
547  {
548  const fieldDescription fieldDesc(is);
549 
550  bool ok = setFaceField::apply(fieldDesc, mesh_, selected_, is);
551 
552  if (!ok)
553  {
554  ok = consumeUnused(fieldDesc, is);
555 
556  if (ok)
557  {
558  // Not meant for us
559  Info<< "Skip " << fieldDesc.type()
560  << " for finite-volume" << nl;
561  }
562  else
563  {
565  << "Unsupported field type: "
566  << fieldDesc.type() << endl;
567  }
568  }
569 
570  return nullptr; // Irrelevant return value
571  }
572  };
573 };
574 
575 
576 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
577 
578 // Dispatcher for setting area face fields
579 struct setAreaField
580 {
581  autoPtr<setAreaField> clone() const { return nullptr; } // placeholder
582 
583  static bool apply
584  (
585  const fieldDescription& fieldDesc,
586  const faMesh& m,
587  const labelList& selectedFaces,
588  Istream& is
589  )
590  {
591  return
592  (
593  setAreaFieldType<scalar>(fieldDesc, m, selectedFaces, is)
594  || setAreaFieldType<vector>(fieldDesc, m, selectedFaces, is)
595  || setAreaFieldType<sphericalTensor>(fieldDesc, m, selectedFaces, is)
596  || setAreaFieldType<symmTensor>(fieldDesc, m, selectedFaces, is)
597  || setAreaFieldType<tensor>(fieldDesc, m, selectedFaces, is)
598  );
599  }
600 
601  class iNew
602  {
603  const faMesh& mesh_;
604  const labelList& selected_;
605 
606  public:
607 
608  iNew(const faMesh& mesh, const labelList& selectedFaces)
609  :
610  mesh_(mesh),
611  selected_(selectedFaces)
612  {}
613 
614  autoPtr<setAreaField> operator()(Istream& is) const
615  {
616  const fieldDescription fieldDesc(is);
617 
618  bool ok = setAreaField::apply(fieldDesc, mesh_, selected_, is);
619 
620  if (!ok)
621  {
622  ok = consumeUnused(fieldDesc, is);
623 
624  if (ok)
625  {
626  // Not meant for us
627  Info<< "Skip " << fieldDesc.type()
628  << " for finite-volume" << nl;
629  }
630  else
631  {
633  << "Unsupported field type: "
634  << fieldDesc.type() << endl;
635  }
636  }
637 
638  return nullptr; // Irrelevant return value
639  }
640  };
641 };
642 
643 
644 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
645 
646 int main(int argc, char *argv[])
647 {
649  (
650  "Set values on a selected set of cells/patch-faces via a dictionary"
651  );
652 
653  argList::addOption("dict", "file", "Alternative setFieldsDict");
654 
656  (
657  "no-finite-area",
658  "Suppress handling of finite-area mesh/fields"
659  );
660 
661  #include "addRegionOption.H"
662  #include "setRootCase.H"
663  #include "createTime.H"
664  #include "createNamedMesh.H"
665 
666  autoPtr<faMesh> faMeshPtr;
667 
668  if (!args.found("no-finite-area"))
669  {
670  faMeshPtr = faMesh::TryNew(mesh);
671  }
672  if (faMeshPtr)
673  {
674  Info<< "Detected finite-area mesh" << nl;
675  }
676 
677  const word dictName("setFieldsDict");
678  #include "setSystemMeshDictionaryIO.H"
679 
680  Info<< "Reading " << dictIO.name() << nl << endl;
681 
682  IOdictionary setFieldsDict(dictIO);
683 
684 
685  // Note.
686  // The PtrList + iNew mechanism is used to trigger the callbacks
687  // and perform the desired actions. The contents of the PtrList
688  // itself are actually irrelevant.
689 
690 
691  // Default field values
692  {
693  const entry* eptr =
694  setFieldsDict.findEntry("defaultFieldValues", keyType::LITERAL);
695 
696  if (eptr)
697  {
698  ITstream& is = eptr->stream();
699 
700  Info<< "Setting volume field default values" << endl;
701 
702  PtrList<setCellField> defaultFieldValues
703  (
704  is,
705  setCellField::iNew(mesh, labelList::null())
706  );
707 
708  if (faMeshPtr)
709  {
710  const faMesh& areaMesh = faMeshPtr();
711  is.rewind();
712 
713  Info<< "Setting area field default values" << endl;
714 
715  PtrList<setAreaField> defaultFieldValues
716  (
717  is,
718  setAreaField::iNew(areaMesh, labelList::null())
719  );
720  }
721  Info<< endl;
722  }
723  }
724 
725 
726  Info<< "Setting field region values" << nl << endl;
727 
728  PtrList<entry> regions(setFieldsDict.lookup("regions"));
729 
730  for (const entry& region : regions)
731  {
732  autoPtr<topoSetSource> source =
733  topoSetSource::New(region.keyword(), mesh, region.dict());
734 
735  if (source().setType() == topoSetSource::CELLSET_SOURCE)
736  {
738  (
739  mesh,
740  "cellSet",
741  mesh.nCells()/10+1 // Reasonable size estimate.
742  );
743 
745 
746  labelList selectedCells(subset.sortedToc());
747 
748  Info<< " Selected "
749  << returnReduce(selectedCells.size(), sumOp<label>())
750  << '/'
752  << " cells" << nl;
753 
754  ITstream& is = region.dict().lookup("fieldValues");
755 
756  PtrList<setCellField> fieldValues
757  (
758  is,
759  setCellField::iNew(mesh, selectedCells)
760  );
761  }
762  else if (source().setType() == topoSetSource::FACESET_SOURCE)
763  {
765  (
766  mesh,
767  "faceSet",
768  mesh.nBoundaryFaces()/10+1
769  );
770 
772 
773  labelList selectedFaces(subset.sortedToc());
774 
775  Info<< " Selected " << selectedFaces.size()
776  << " faces" << nl;
777 
778  ITstream& is = region.dict().lookup("fieldValues");
779 
780  PtrList<setFaceField> fieldValues
781  (
782  is,
783  setFaceField::iNew(mesh, selectedFaces)
784  );
785 
786  if (faMeshPtr)
787  {
788  const faMesh& areaMesh = faMeshPtr();
789 
790  const labelUList& faceLabels = areaMesh.faceLabels();
791 
792  // Transcribe from mesh faces to finite-area addressing
793  labelList areaFaces(faceLabels.size());
794 
795  label nUsed = 0;
796  forAll(faceLabels, facei)
797  {
798  const label meshFacei = faceLabels[facei];
799 
800  if (subset.test(meshFacei))
801  {
802  areaFaces[nUsed] = facei;
803  ++nUsed;
804  }
805  }
806  areaFaces.resize(nUsed);
807 
808  Info<< " Selected "
809  << returnReduce(areaFaces.size(), sumOp<label>())
810  << '/'
811  << returnReduce(faceLabels.size(), sumOp<label>())
812  << " area faces" << nl;
813 
814  is.rewind();
815 
816  PtrList<setAreaField> fieldValues
817  (
818  is,
819  setAreaField::iNew(areaMesh, areaFaces)
820  );
821  }
822  }
823  }
824 
825 
826  Info<< "\nEnd\n" << endl;
827 
828  return 0;
829 }
830 
831 
832 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:88
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
rDeltaTY field()
A list of face labels.
Definition: faceSet.H:47
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
Create a new set and ADD elements to it.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
const word dictName("faMeshDefinition")
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
const labelList & patchID() const
Per boundary face label the patch index.
Required Variables.
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:365
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
static const List< label > & null()
Return a null List.
Definition: ListI.H:102
SubField is a Field obtained as a section of another Field, without its own allocation. SubField is derived from a SubList rather than a List.
Definition: Field.H:62
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:377
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
Generic templated field type.
Definition: Field.H:61
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Time & time() const noexcept
Return time registry.
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:376
label nInternalFaces() const noexcept
Number of internal faces.
String literal.
Definition: keyType.H:82
virtual void rewind()
Rewind the stream so that it may be read again.
Definition: ITstream.C:525
const direction noexcept
Definition: Scalar.H:258
bool isNull(const T *ptr)
True if ptr is a pointer (of type T) to the nullObject.
Definition: nullObject.H:227
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:89
static autoPtr< topoSetSource > New(const word &topoSetSourceType, const polyMesh &mesh, const dictionary &dict)
Return a reference to the selected topoSetSource.
messageStream Warning
Warning stream (stdout output on master, null elsewhere), with additional &#39;FOAM Warning&#39; header text...
List< T > subset(const BoolListType &select, const UList< T > &input, const bool invert=false)
Extract elements of the input list when select is true.
virtual void applyToSet(const topoSetSource::setAction action, topoSet &set) const =0
Apply specified action to the topoSet.
#define WarningInFunction
Report a warning using Foam::Warning.
label nCells() const noexcept
Number of mesh cells.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
A collection of cell labels.
Definition: cellSet.H:47
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
static autoPtr< faMesh > TryNew(const polyMesh &pMesh)
Read construction from polyMesh if all files are available.
Definition: faMeshNew.C:149
Mesh data needed to do the Finite Area discretisation.
Definition: areaFaMesh.H:47
messageStream Info
Information stream (stdout output on master, null elsewhere)
IOobject dictIO
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual ITstream & stream() const =0
Return token stream, if entry is a primitive entry.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
bool found
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
An input stream of tokens.
Definition: ITstream.H:48
Namespace for OpenFOAM.
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
After completion all processors have the same data.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157