sampledSets.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 OpenFOAM Foundation
9  Copyright (C) 2015-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 \*---------------------------------------------------------------------------*/
28 
29 #include "sampledSets.H"
30 #include "dictionary.H"
31 #include "Time.H"
32 #include "globalIndex.H"
33 #include "volFields.H"
34 #include "mapPolyMesh.H"
35 #include "IOmanip.H"
36 #include "IOobjectList.H"
37 #include "IndirectList.H"
38 #include "ListOps.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(sampledSets, 0);
46 
48  (
49  functionObject,
50  sampledSets,
51  dictionary
52  );
53 }
54 
55 #include "sampledSetsImpl.C"
56 
57 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58 
59 Foam::autoPtr<Foam::coordSetWriter> Foam::sampledSets::newWriter
60 (
61  word writerType,
62  const dictionary& topDict,
63  const dictionary& setDict
64 )
65 {
66  // Per-set adjustment
67  setDict.readIfPresent<word>("setFormat", writerType);
68 
69  return coordSetWriter::New
70  (
71  writerType,
72  // Top-level/set-specific "formatOptions"
73  coordSetWriter::formatOptions(topDict, setDict, writerType)
74  );
75 }
76 
77 
78 Foam::OFstream* Foam::sampledSets::createProbeFile(const word& fieldName)
79 {
80  // Open new output stream
81 
82  OFstream* osptr = probeFilePtrs_.lookup(fieldName, nullptr);
83 
84  if (!osptr && Pstream::master())
85  {
86  // Put in undecomposed case
87  // (Note: gives problems for distributed data running)
88 
89  fileName probeDir
90  (
91  mesh_.time().globalPath()
93  / name()/mesh_.regionName()
94  // Use startTime as the instance for output files
95  / mesh_.time().timeName(mesh_.time().startTime().value())
96  );
97  probeDir.clean(); // Remove unneeded ".."
98 
99  // Create directory if needed
100  Foam::mkDir(probeDir);
101 
102  probeFilePtrs_.insert
103  (
104  fieldName,
105  autoPtr<OFstream>::New(probeDir/fieldName)
106  );
107  osptr = probeFilePtrs_.lookup(fieldName, nullptr);
108 
109  if (osptr)
110  {
111  auto& os = *osptr;
112 
113  DebugInfo<< "open probe stream: " << os.name() << endl;
114 
115  const unsigned int width(IOstream::defaultPrecision() + 7);
116 
117  label nPoints = 0;
118  forAll(*this, seti)
119  {
120  const coordSet& s = gatheredSets_[seti];
121 
122  const pointField& pts = static_cast<const pointField&>(s);
123 
124  for (const point& p : pts)
125  {
126  os << "# Probe " << nPoints++ << ' ' << p << nl;
127  }
128  }
129 
130  os << '#' << setw(IOstream::defaultPrecision() + 6)
131  << "Probe";
132 
133  for (label probei = 0; probei < nPoints; ++probei)
134  {
135  os << ' ' << setw(width) << probei;
136  }
137  os << nl;
138 
139  os << '#' << setw(IOstream::defaultPrecision() + 6)
140  << "Time" << endl;
141  }
142  }
143 
144  return osptr;
145 }
146 
147 
148 void Foam::sampledSets::gatherAllSets()
149 {
150  // Any writer references will become invalid
151  for (auto& writer : writers_)
152  {
153  writer.expire();
154  }
155 
156  const PtrList<sampledSet>& localSets = *this;
157 
158  gatheredSets_.clear();
159  gatheredSets_.resize(localSets.size());
160  gatheredSorting_.resize_nocopy(localSets.size());
161  globalIndices_.resize_nocopy(localSets.size());
162 
163  forAll(localSets, seti)
164  {
165  const coordSet& coords = localSets[seti];
166 
167  globalIndices_[seti].reset(globalIndex::gatherOnly{}, coords.size());
168  gatheredSets_.set(seti, coords.gatherSort(gatheredSorting_[seti]));
169  }
170 }
171 
172 
173 Foam::IOobjectList Foam::sampledSets::preCheckFields(unsigned request)
174 {
175  wordList allFields; // Just needed for warnings
176  HashTable<wordHashSet> selected;
177 
178  IOobjectList objects(0);
179 
180  if (loadFromFiles_)
181  {
182  // Check files for a particular time
183  objects = IOobjectList(mesh_, mesh_.time().timeName());
184 
185  allFields = objects.names();
186  selected = objects.classes(fieldSelection_);
187  }
188  else
189  {
190  // Check currently available fields
191  allFields = mesh_.names();
192  selected = mesh_.classes(fieldSelection_);
193  }
194 
195  // Parallel consistency (no-op in serial)
196  // Probably not needed...
198 
199 
200  DynamicList<label> missed(fieldSelection_.size());
201 
202  // Detect missing fields
203  forAll(fieldSelection_, i)
204  {
205  if
206  (
207  fieldSelection_[i].isLiteral()
208  && !ListOps::found(allFields, fieldSelection_[i])
209  )
210  {
211  missed.append(i);
212  }
213  }
214 
215  if (missed.size() && (request != ACTION_NONE))
216  {
218  << nl
219  << "Cannot find "
220  << (loadFromFiles_ ? "field file" : "registered field")
221  << " matching "
222  << UIndirectList<wordRe>(fieldSelection_, missed) << endl;
223  }
224 
225 
226  // The selected field names, ordered by (scalar, vector, ...)
227  // with internal sorting
228 
229  selectedFieldNames_.clear();
230 
231  do
232  {
233  #undef doLocalCode
234  #define doLocalCode(InputType) \
235  { \
236  const auto iter = selected.find(InputType::typeName); \
237  if (iter.found()) \
238  { \
239  selectedFieldNames_.append(iter.val().sortedToc()); \
240  } \
241  }
242 
248  #undef doLocalCode
249  }
250  while (false);
251 
252 
253  // Now propagate field counts (per surface)
254  // - can update writer even when not writing without problem
255 
256  const label nFields = selectedFieldNames_.size();
257 
258  if (writeAsProbes_)
259  {
260  // Close streams for fields that no longer exist
261  forAllIters(probeFilePtrs_, iter)
262  {
263  if (!selectedFieldNames_.found(iter.key()))
264  {
265  DebugInfo
266  << "close probe stream: "
267  << iter()->name() << endl;
268 
269  probeFilePtrs_.remove(iter);
270  }
271  }
272  }
273  else if ((request & ACTION_WRITE) != 0)
274  {
275  forAll(writers_, seti)
276  {
277  coordSetWriter& writer = writers_[seti];
278 
279  writer.nFields(nFields);
280  }
281  }
282 
283  return objects;
284 }
285 
286 
287 void Foam::sampledSets::initDict(const dictionary& dict, const bool initial)
288 {
290  if (initial)
291  {
292  writers_.clear();
293  }
294 
295  const entry* eptr = dict.findEntry("sets");
296 
297  if (eptr && eptr->isDict())
298  {
299  PtrList<sampledSet> sampSets(eptr->dict().size());
300  if (initial && !writeAsProbes_)
301  {
302  writers_.resize(sampSets.size());
303  }
304 
305  label seti = 0;
306 
307  for (const entry& dEntry : eptr->dict())
308  {
309  if (!dEntry.isDict())
310  {
311  continue;
312  }
313 
314  const dictionary& subDict = dEntry.dict();
315 
316  autoPtr<sampledSet> sampSet =
318  (
319  dEntry.keyword(),
320  mesh_,
321  searchEngine_,
322  subDict
323  );
324 
325  // if (!sampSet || !sampSet->enabled())
326  // {
327  // continue;
328  // }
329 
330  // Define the set
331  sampSets.set(seti, sampSet);
332 
333  // Define writer, but do not attached
334  if (initial && !writeAsProbes_)
335  {
336  writers_.set
337  (
338  seti,
339  newWriter(writeFormat_, dict_, subDict)
340  );
341 
342  // Use outputDir/TIME/set-name
343  writers_[seti].useTimeDir(true);
344  writers_[seti].verbose(verbose_);
345  }
346  ++seti;
347  }
348 
349  sampSets.resize(seti);
350  if (initial && !writeAsProbes_)
351  {
352  writers_.resize(seti);
353  }
354  static_cast<PtrList<sampledSet>&>(*this).transfer(sampSets);
355  }
356  else if (eptr)
357  {
358  // This is slightly trickier.
359  // We want access to the individual dictionaries used for construction
360 
361  DynamicList<dictionary> capture;
362 
363  PtrList<sampledSet> input
364  (
365  eptr->stream(),
366  sampledSet::iNewCapture(mesh_, searchEngine_, capture)
367  );
368 
369  PtrList<sampledSet> sampSets(input.size());
370  if (initial && !writeAsProbes_)
371  {
372  writers_.resize(sampSets.size());
373  }
374 
375  label seti = 0;
376 
377  forAll(input, inputi)
378  {
379  const dictionary& subDict = capture[inputi];
380 
381  autoPtr<sampledSet> sampSet = input.release(inputi);
382 
383  // if (!sampSet || !sampSet->enabled())
384  // {
385  // continue;
386  // }
387 
388  // Define the set
389  sampSets.set(seti, sampSet);
390 
391  // Define writer, but do not attached
392  if (initial && !writeAsProbes_)
393  {
394  writers_.set
395  (
396  seti,
397  newWriter(writeFormat_, dict_, subDict)
398  );
399 
400  // Use outputDir/TIME/set-name
401  writers_[seti].useTimeDir(true);
402  writers_[seti].verbose(verbose_);
403  }
404  ++seti;
405  }
406 
407  sampSets.resize(seti);
408  if (initial && !writeAsProbes_)
409  {
410  writers_.resize(seti);
411  }
412 
413  static_cast<PtrList<sampledSet>&>(*this).transfer(sampSets);
414  }
415 
416  gatherAllSets();
417 
418  needsCorrect_ = false;
419 }
420 
421 
422 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
423 
424 Foam::sampledSets::sampledSets
425 (
426  const word& name,
427  const Time& runTime,
428  const dictionary& dict
429 )
430 :
431  functionObjects::fvMeshFunctionObject(name, runTime, dict),
432  PtrList<sampledSet>(),
433  dict_(dict),
434  loadFromFiles_(false),
435  verbose_(false),
436  onExecute_(false),
437  needsCorrect_(false),
438  writeAsProbes_(false),
439  outputPath_
440  (
441  time_.globalPath()/functionObject::outputPrefix
442  / name/mesh_.regionName()
443  ),
444  searchEngine_(mesh_),
445  samplePointScheme_(),
446  writeFormat_(),
447  selectedFieldNames_(),
448  writers_(),
449  probeFilePtrs_(),
450  gatheredSets_(),
451  gatheredSorting_(),
452  globalIndices_()
453 {
454  outputPath_.clean(); // Remove unneeded ".."
455 
456  read(dict);
457 }
458 
459 
460 Foam::sampledSets::sampledSets
461 (
462  const word& name,
463  const objectRegistry& obr,
464  const dictionary& dict,
465  const bool loadFromFiles
466 )
467 :
468  functionObjects::fvMeshFunctionObject(name, obr, dict),
469  PtrList<sampledSet>(),
470  dict_(dict),
471  loadFromFiles_(loadFromFiles),
472  verbose_(false),
473  onExecute_(false),
474  needsCorrect_(false),
475  writeAsProbes_(false),
476  outputPath_
477  (
478  time_.globalPath()/functionObject::outputPrefix
479  / name/mesh_.regionName()
480  ),
481  searchEngine_(mesh_),
482  samplePointScheme_(),
483  writeFormat_(),
484  selectedFieldNames_(),
485  writers_(),
486  probeFilePtrs_(),
487  gatheredSets_(),
488  gatheredSorting_(),
489  globalIndices_()
490 {
491  outputPath_.clean(); // Remove unneeded ".."
493  read(dict);
494 }
495 
496 
497 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
498 
499 bool Foam::sampledSets::verbose(const bool on) noexcept
500 {
501  bool old(verbose_);
502  verbose_ = on;
503  return old;
504 }
505 
506 
508 {
509  if (&dict_ != &dict)
510  {
511  // Update local copy of dictionary
512  dict_ = dict;
513  }
514 
516 
518  writers_.clear();
519  fieldSelection_.clear();
520  selectedFieldNames_.clear();
521 
522  gatheredSets_.clear();
523  gatheredSorting_.clear();
524  globalIndices_.clear();
525 
526  verbose_ = dict.getOrDefault("verbose", false);
527  onExecute_ = dict.getOrDefault("sampleOnExecute", false);
528 
529  samplePointScheme_ =
530  dict.getOrDefault<word>("interpolationScheme", "cellPoint");
531 
532  const entry* eptr = dict.findEntry("sets");
533 
534  if (eptr)
535  {
536  dict.readEntry("setFormat", writeFormat_);
537  }
538 
539  // Hard-coded handling of ensemble 'probes' writer
540  writeAsProbes_ = ("probes" == writeFormat_);
541  if (!writeAsProbes_)
542  {
543  // Close all streams
544  probeFilePtrs_.clear();
545  }
546 
547  initDict(dict, true);
548 
549  // Have some sets, so sort out which fields are needed and report
550 
551  if (this->size())
552  {
553  dict_.readEntry("fields", fieldSelection_);
554  fieldSelection_.uniq();
555 
556  // Report
557  if (writeAsProbes_)
558  {
559  Info<< "Sampled set as probes ensemble:" << nl;
560 
561  forAll(*this, seti)
562  {
563  const sampledSet& s = (*this)[seti];
564  Info<< " " << s.name();
565  }
566  Info<< nl;
567  }
568  else
569  {
570  Info<< "Sampled set:" << nl;
571 
572  forAll(*this, seti)
573  {
574  const sampledSet& s = (*this)[seti];
575 
576  Info<< " " << s.name() << " -> "
577  << writers_[seti].type() << nl;
578  }
579  }
580 
581  Info<< endl;
582  }
583 
584  if (debug && Pstream::master())
585  {
586  Pout<< "sample fields:" << flatOutput(fieldSelection_) << nl
587  << "sample sets:" << nl << '(' << nl;
588 
589  for
590  (
591  const sampledSet& s
592  : static_cast<const PtrList<sampledSet>&>(*this)
593  )
594  {
595  Pout<< " " << s << endl;
596  }
597  Pout<< ')' << endl;
598  }
599 
600  if (writeAsProbes_)
601  {
602  (void) preCheckFields(ACTION_NONE);
603  }
604 
605  // FUTURE:
606  // Ensure all sets and merge information are expired
607  // expire(true);
608 
609  return true;
610 }
611 
612 
613 bool Foam::sampledSets::performAction(unsigned request)
614 {
615  if (empty())
616  {
617  // Nothing to do
618  return true;
619  }
620  else if (needsCorrect_)
621  {
622  searchEngine_.correct();
623  initDict(dict_, false);
624  }
625 
626  // FUTURE:
627  // Update sets and store
628  // ...
629 
630  // Determine availability of fields.
631  // Count number of fields (only seems to be needed for VTK legacy)
632 
633  IOobjectList objects = preCheckFields(request);
634 
635  const label nFields = selectedFieldNames_.size();
636 
637  if (!nFields)
638  {
639  // Nothing to do
640  return true;
641  }
642 
643  // Update writers
644  if (!writeAsProbes_)
645  {
646  forAll(*this, seti)
647  {
648  const coordSet& s = gatheredSets_[seti];
649 
650  if ((request & ACTION_WRITE) != 0)
651  {
652  coordSetWriter& writer = writers_[seti];
653 
654  if (writer.needsUpdate())
655  {
656  writer.setCoordinates(s);
657  }
658 
659  if (writer.buffering())
660  {
661  writer.open
662  (
663  outputPath_
664  / word
665  (
666  s.name()
667  + coordSetWriter::suffix(selectedFieldNames_)
668  )
669  );
670  }
671  else
672  {
673  writer.open(outputPath_/s.name());
674  }
675 
676  writer.beginTime(mesh_.time());
677  }
678  }
679  }
680 
681  // Sample fields
682 
683  performAction<VolumeField<scalar>>(objects, request);
684  performAction<VolumeField<vector>>(objects, request);
685  performAction<VolumeField<sphericalTensor>>(objects, request);
686  performAction<VolumeField<symmTensor>>(objects, request);
687  performAction<VolumeField<tensor>>(objects, request);
688 
689 
690  // Finish this time step
691  if (!writeAsProbes_)
692  {
693  forAll(writers_, seti)
694  {
695  // Write geometry if no fields were written so that we still
696  // can have something to look at
697 
698  if ((request & ACTION_WRITE) != 0)
699  {
704 
705  writers_[seti].endTime();
706  }
707  }
708  }
709 
710  return true;
711 }
712 
713 
715 {
716  if (onExecute_)
717  {
718  return performAction(ACTION_ALL & ~ACTION_WRITE);
719  }
720 
721  return true;
722 }
723 
726 {
727  return performAction(ACTION_ALL);
728 }
729 
732 {
733  needsCorrect_ = true;
734 }
735 
736 
738 {
739  if (&mpm.mesh() == &mesh_)
740  {
741  correct();
742  }
743 }
744 
745 
746 void Foam::sampledSets::movePoints(const polyMesh& mesh)
747 {
748  if (&mesh == &mesh_)
749  {
750  correct();
751  }
752 }
753 
754 
756 {
757  if (state != polyMesh::UNCHANGED)
758  {
759  correct();
760  }
761 }
762 
763 
764 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:88
dictionary dict
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: sampledSets.C:739
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
static autoPtr< sampledSet > New(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Return a reference to the selected sampledSet.
Definition: sampledSet.C:477
bool found(const ListType &input, const UnaryPredicate &pred, const label start=0)
Same as found_if.
Definition: ListOps.H:821
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:89
Output to file stream, using an OSstream.
Definition: OFstream.H:49
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary.
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:416
virtual bool read(const dictionary &)
Read the sampledSets.
Definition: sampledSets.C:500
Abstract base-class for Time/database function objects.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
Definition: volFieldsFwd.H:87
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:157
Macros for easy insertion into run-time selection tables.
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
static word suffix(const word &fldName, const word &fileExt=word::null)
Name suffix based on fieldName (underscore separator)
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
static autoPtr< coordSetWriter > New(const word &writeFormat)
Return a reference to the selected writer.
Foam::word regionName(Foam::polyMesh::defaultRegion)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
#define doLocalCode(InputType)
dynamicFvMesh & mesh
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:567
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:79
A class for handling words, derived from Foam::string.
Definition: word.H:63
label nPoints
static Istream & input(Istream &is, IntRange< T > &range)
Definition: IntRanges.C:48
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:328
void correct()
Correct for mesh changes.
Definition: sampledSets.C:724
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: sampledSets.C:730
bool verbose(const bool on) noexcept
Enable/disable verbose output.
Definition: sampledSets.C:492
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
#define DebugInfo
Report an information message using Foam::Info.
virtual bool execute()
Execute, currently does nothing.
Definition: sampledSets.C:707
const direction noexcept
Definition: Scalar.H:258
Istream and Ostream manipulators taking arguments.
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
int debug
Static debugging option.
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2212/OpenFOAM-v2212/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
OBJstream os(runTime.globalPath()/outputName)
virtual bool open(const fileName &file, bool parallel=UPstream::parRun())
Open file for writing (creates parent directory).
defineTypeNameAndDebug(combustionModel, 0)
static bool clean(std::string &str)
Cleanup filename string, possibly applies other transformations such as changing the path separator e...
Definition: fileName.C:192
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:128
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:58
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
Definition: sampledSets.C:748
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
virtual bool write()
Sample and write.
Definition: sampledSets.C:718
const polyMesh & mesh() const
Return polyMesh.
Definition: mapPolyMesh.H:437
static bool master(const label communicator=worldComm)
Am I the master rank.
Definition: UPstream.H:672
static word outputPrefix
Directory prefix.
const entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:80
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:90
messageStream Info
Information stream (stdout output on master, null elsewhere)
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:89
volScalarField & p
Registry of regIOobjects.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
const pointField & pts