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