surfaceWriter.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) 2019-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "surfaceWriter.H"
29 #include "proxySurfaceWriter.H"
30 #include "MeshedSurfaceProxy.H"
31 
32 #include "fileFormats.H"
33 #include "Time.H"
34 #include "coordinateRotation.H"
35 #include "transformField.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(surfaceWriter, 0);
43  defineRunTimeSelectionTable(surfaceWriter, word);
44  defineRunTimeSelectionTable(surfaceWriter, wordDict);
45 }
46 
47 Foam::scalar Foam::surfaceWriter::defaultMergeDim = 1e-8;
48 
49 
50 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
51 
52 bool Foam::surfaceWriter::supportedType(const word& writeType)
53 {
54  return
55  (
56  wordConstructorTablePtr_->found(writeType)
57  || wordDictConstructorTablePtr_->found(writeType)
59  );
60 }
61 
62 
64 (
65  const dictionary& dict,
66  const word& formatName,
67  const word& entryName
68 )
69 {
70  return fileFormats::getFormatOptions(dict, formatName, entryName);
71 }
72 
73 
75 (
76  const dictionary& dict,
77  const dictionary& surfDict,
78  const word& formatName,
79  const word& entryName
80 )
81 {
82  return fileFormats::getFormatOptions(dict, surfDict, formatName, entryName);
83 }
84 
85 
87 Foam::surfaceWriter::New(const word& writeType)
88 {
89  // Constructors without dictionary options
90  auto* ctorPtr = wordConstructorTable(writeType);
91 
92  if (!ctorPtr)
93  {
95  {
96  // Generally unknown, but handle via 'proxy' handler
98  (
99  new surfaceWriters::proxyWriter(writeType)
100  );
101  }
102 
104  << "Unknown write type \"" << writeType << "\"\n\n"
105  << "Valid write types : "
106  << flatOutput(wordConstructorTablePtr_->sortedToc()) << nl
107  << "Valid proxy types : "
109  << exit(FatalError);
110  }
112  return autoPtr<surfaceWriter>(ctorPtr());
113 }
114 
115 
118 (
119  const word& writeType,
120  const dictionary& writeOpts
121 )
122 {
123  // Constructors with dictionary options
124  {
125  auto* ctorPtr = wordDictConstructorTable(writeType);
126 
127  if (ctorPtr)
128  {
129  return autoPtr<surfaceWriter>(ctorPtr(writeOpts));
130  }
131  }
132 
133 
134  // Constructors without dictionary options
135  auto* ctorPtr = wordConstructorTable(writeType);
136 
137  if (!ctorPtr)
138  {
140  {
141  // Generally unknown, but handle via 'proxy' handler
142  return autoPtr<surfaceWriter>
143  (
144  new surfaceWriters::proxyWriter(writeType, writeOpts)
145  );
146  }
147 
149  << "Unknown write type \"" << writeType << "\"\n\n"
150  << "Valid write types : "
151  << wordConstructorTablePtr_->sortedToc() << nl
152  << "Valid proxy types : "
154  << exit(FatalError);
155  }
157  return autoPtr<surfaceWriter>(ctorPtr());
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
162 
164 :
165  surf_(),
166  mergedSurf_(),
167  adjustedSurf_(),
168  mergeDim_(defaultMergeDim),
169  geometryScale_(1),
170  geometryCentre_(Zero),
171  geometryTransform_(),
172  upToDate_(false),
173  wroteGeom_(false),
174  parallel_(true),
175  useTimeDir_(false),
176  isPointData_(false),
177  verbose_(false),
178  commType_(UPstream::commsTypes::scheduled),
179  nFields_(0),
180  currTime_(),
181  outputPath_(),
182  fieldLevel_(),
183  fieldScale_()
184 {
186 }
187 
188 
190 :
191  surfaceWriter()
192 {
193  options.readIfPresent("verbose", verbose_);
194 
195  UPstream::commsTypeNames.readIfPresent("commsType", options, commType_);
196 
197  geometryScale_ = 1;
200 
201  options.readIfPresent("scale", geometryScale_);
202 
203  // Optional cartesian coordinate system transform
204  const auto* dictptr = options.findDict("transform", keyType::LITERAL);
205 
206  if (dictptr)
207  {
208  dictptr->readIfPresent("rotationCentre", geometryCentre_);
209 
210  // 'origin' is optional within sub-dictionary
213  }
214 
215  fieldLevel_ = options.subOrEmptyDict("fieldLevel");
216  fieldScale_ = options.subOrEmptyDict("fieldScale");
217 
218  if (verbose_)
219  {
220  Info<< "Create surfaceWriter ("
221  << (this->isPointData() ? "point" : "face") << " data):"
222  << " commsType=" << UPstream::commsTypeNames[commType_] << endl;
223  }
224 }
225 
226 
227 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
228 
230 {
231  close();
232 }
233 
234 
235 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
237 void Foam::surfaceWriter::setTime(const instant& inst)
238 {
239  currTime_ = inst;
240 }
241 
243 void Foam::surfaceWriter::setTime(scalar timeValue)
244 {
245  currTime_ = instant(timeValue);
246 }
247 
248 
249 void Foam::surfaceWriter::setTime(scalar timeValue, const word& timeName)
250 {
251  currTime_.value() = timeValue;
252  currTime_.name() = timeName;
253 }
254 
255 
257 {
258  currTime_.value() = 0;
259  currTime_.name().clear();
260 }
261 
264 {
265  setTime(t.value(), t.timeName());
266 }
267 
269 void Foam::surfaceWriter::beginTime(const instant& inst)
270 {
271  setTime(inst);
272 }
273 
276 {
277  unsetTime();
278 }
279 
280 
281 void Foam::surfaceWriter::open(const fileName& outputPath)
282 {
283  outputPath_ = outputPath;
284  wroteGeom_ = false;
285 }
286 
287 
289 (
290  const meshedSurf& surf,
291  const fileName& outputPath,
292  bool parallel
293 )
294 {
295  close();
296  setSurface(surf, parallel);
297  open(outputPath);
298 }
299 
300 
302 (
303  const pointField& points,
304  const faceList& faces,
305  const fileName& outputPath,
306  bool parallel
307 )
308 {
309  close();
310  setSurface(points, faces, parallel);
311  open(outputPath);
312 }
313 
314 
316 (
317  const meshedSurf& surf,
318  const fileName& outputPath
319 )
320 {
321  close();
322  setSurface(surf, parallel_);
323  open(outputPath);
324 }
325 
326 
328 (
329  const pointField& points,
330  const faceList& faces,
331  const fileName& outputPath
332 )
333 {
334  close();
335  setSurface(points, faces, parallel_);
336  open(outputPath);
337 }
338 
339 
341 {
342  outputPath_.clear();
343  wroteGeom_ = false;
344 }
345 
346 
348 {
349  close();
350  expire();
351  surf_.clear();
352 }
353 
354 
356 (
357  const meshedSurf& surf,
358  bool parallel
359 )
360 {
361  expire();
362  surf_.reset(surf);
363  parallel_ = (parallel && Pstream::parRun());
364 }
365 
366 
368 (
369  const pointField& points,
370  const faceList& faces,
371  bool parallel
372 )
373 {
374  expire();
375  surf_.reset(points, faces);
376  parallel_ = (parallel && Pstream::parRun());
377 }
378 
379 
381 (
382  const meshedSurf& surf
383 )
384 {
385  setSurface(surf, parallel_);
386 }
387 
388 
390 (
391  const pointField& points,
392  const faceList& faces
393 )
394 {
395  setSurface(points, faces, parallel_);
396 }
397 
400 {
401  return !upToDate_;
402 }
403 
406 {
407  return wroteGeom_;
408 }
409 
410 
412 {
413  const bool changed = upToDate_;
414 
415  upToDate_ = false;
416  wroteGeom_ = false;
417  adjustedSurf_.clear();
418  mergedSurf_.clear();
419 
420  // Field count (nFields_) is a different type of accounting
421  // and is unaffected by geometry changes
422 
423  return changed;
424 }
425 
428 {
429  return surf_.good();
430 }
431 
432 
433 bool Foam::surfaceWriter::empty() const
434 {
435  const bool value = surf_.faces().empty();
436 
437  return (parallel_ ? returnReduceAnd(value) : value);
438 }
439 
440 
441 Foam::label Foam::surfaceWriter::size() const
442 {
443  const label value = surf_.faces().size();
444 
445  return (parallel_ ? returnReduce(value, sumOp<label>()) : value);
446 }
447 
448 
450 {
451  if (!is_open())
452  {
454  << type() << " : Attempted to write without a path" << nl
455  << exit(FatalError);
456  }
457 }
458 
459 
460 bool Foam::surfaceWriter::merge() const
461 {
462  bool changed = false;
463 
464  if (!upToDate_)
465  {
466  // Similar to expire
467  adjustedSurf_.clear();
468 
469  if (parallel_ && Pstream::parRun())
470  {
471  changed = mergedSurf_.merge(surf_, mergeDim_);
472  }
473  else
474  {
475  mergedSurf_.clear();
476  }
477  }
478 
479  upToDate_ = true;
480 
481  if (changed)
482  {
483  wroteGeom_ = false;
484  }
485 
486  return changed;
487 }
488 
489 
491 {
492  merge();
493 
494  if (parallel_ && Pstream::parRun())
495  {
496  return mergedSurf_;
497  }
498 
499  return surf_;
500 }
501 
502 
504 {
505  if (!upToDate_)
506  {
507  adjustedSurf_.clear();
508  }
509 
510  if (!adjustedSurf_.good())
511  {
512  adjustedSurf_.reset(surface());
513 
514  tmp<pointField> tpts;
515 
516  if (geometryTransform_.valid())
517  {
518  if (!geometryTransform_.R().is_identity())
519  {
520  if (magSqr(geometryCentre_) > ROOTVSMALL)
521  {
522  // Set centre of rotation,
523  // followed by forward transform (local -> global)
524  tpts =
525  geometryTransform_.globalPosition
526  (
527  adjustedSurf_.points0() - geometryCentre_
528  );
529 
530  // Unset centre of rotation
531  tpts.ref() += geometryCentre_;
532  }
533  else
534  {
535  // Forward transform (local -> global)
536  tpts =
537  geometryTransform_.globalPosition
538  (
539  adjustedSurf_.points0()
540  );
541  }
542  }
543  else if (magSqr(geometryTransform_.origin()) > ROOTVSMALL)
544  {
545  // Translate only (local -> global)
546  tpts = (adjustedSurf_.points0() + geometryTransform_.origin());
547  }
548  }
549 
550  adjustedSurf_.movePoints(tpts);
551  adjustedSurf_.scalePoints(geometryScale_);
552  }
553 
554  return adjustedSurf_;
555 }
556 
557 
558 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
559 
560 template<class Type>
562 (
563  const Field<Type>& fld
564 ) const
565 {
566  if (parallel_ && Pstream::parRun())
567  {
568  // Ensure geometry is also merged
569  merge();
570 
571  // Gather all values
572  auto tfield = tmp<Field<Type>>::New();
573  auto& allFld = tfield.ref();
574 
575  // Update any expired global index (as required)
576 
577  const globalIndex& globIndex =
578  (
579  this->isPointData()
580  ? mergedSurf_.pointGlobalIndex()
581  : mergedSurf_.faceGlobalIndex()
582  );
583 
584  globIndex.gather
585  (
586  fld,
587  allFld,
589  commType_,
591  );
592 
593  // Renumber (point data) to correspond to merged points
594  if
595  (
597  && this->isPointData()
598  && mergedSurf_.pointsMap().size()
599  )
600  {
601  inplaceReorder(mergedSurf_.pointsMap(), allFld);
602  allFld.resize(mergedSurf_.points().size());
603  }
604 
605  return tfield;
606  }
607 
608  // Mark that any geometry changes have been taken care of
609  upToDate_ = true;
611  return fld;
612 }
613 
614 
615 template<class Type>
617 (
618  const word& fieldName,
619  const tmp<Field<Type>>& tfield
620 ) const
621 {
622  if (verbose_)
623  {
624  Info<< "Writing field " << fieldName;
625  }
626 
627  tmp<Field<Type>> tadjusted;
628 
629  // Output scaling for the variable, but not for integer types
630  // which are typically ids etc.
631  if (!std::is_integral<Type>::value)
632  {
633  scalar value;
634 
635  // Remove *uniform* reference level
636  if
637  (
638  fieldLevel_.readIfPresent(fieldName, value, keyType::REGEX)
639  && !equal(value, 0)
640  )
641  {
642  // Could also detect brackets (...) and read accordingly
643  // or automatically scale by 1/sqrt(nComponents) instead ...
644 
645  Type refLevel;
646  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; ++cmpt)
647  {
648  setComponent(refLevel, cmpt) = value;
649  }
650 
651  if (verbose_)
652  {
653  Info<< " [level " << refLevel << ']';
654  }
655 
656  if (!tadjusted)
657  {
658  // Steal or clone
659  tadjusted.reset(tfield.ptr());
660  }
661 
662  // Remove offset level
663  tadjusted.ref() -= refLevel;
664  }
665 
666  // Apply scaling
667  if
668  (
669  fieldScale_.readIfPresent(fieldName, value, keyType::REGEX)
670  && !equal(value, 1)
671  )
672  {
673  if (verbose_)
674  {
675  Info<< " [scaling " << value << ']';
676  }
677 
678  if (!tadjusted)
679  {
680  // Steal or clone
681  tadjusted.reset(tfield.ptr());
682  }
683 
684  // Apply scaling
685  tadjusted.ref() *= value;
686  }
687 
688  // Rotate fields (vector and non-spherical tensors)
689  if
690  (
691  (pTraits<Type>::rank != 0 && pTraits<Type>::nComponents > 1)
692  && geometryTransform_.valid()
693  && !geometryTransform_.R().is_identity()
694  )
695  {
696  if (!tadjusted)
697  {
698  // Steal or clone
699  tadjusted.reset(tfield.ptr());
700  }
701 
703  (
704  tadjusted.ref(),
705  geometryTransform_.R(),
706  tadjusted()
707  );
708  }
709  }
710 
711  return (tadjusted ? tadjusted : tfield);
712 }
713 
714 
715 #define defineSurfaceFieldMethods(ThisClass, Type) \
716  Foam::tmp<Foam::Field<Type>> \
717  ThisClass::mergeField(const Field<Type>& fld) const \
718  { \
719  return mergeFieldTemplate(fld); \
720  } \
721  \
722  Foam::tmp<Foam::Field<Type>> \
723  ThisClass::adjustField \
724  ( \
725  const word& fieldName, \
726  const tmp<Field<Type>>& tfield \
727  ) const \
728  { \
729  return adjustFieldTemplate(fieldName, tfield); \
730  }
731 
739 #undef defineSurfaceFieldMethod
740 
741 
742 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
743 
744 Foam::Ostream& Foam::operator<<
745 (
746  Ostream& os,
747  const InfoProxy<surfaceWriter>& iproxy
748 )
749 {
750  const auto& w = *iproxy;
751 
752  os << "surfaceWriter:"
753  << " upToDate: " << w.upToDate_
754  << " PointData: " << w.isPointData_
755  << " nFields: " << w.nFields_
756  << " time: " << w.currTime_
757  << " path: " << w.outputPath_ << endl;
758 
759  return os;
760 }
761 
762 
763 // ************************************************************************* //
const Type & value() const noexcept
Return const reference to value.
dictionary dict
label size() const
The global number of faces for the associated surface.
virtual bool expire()
Mark that surface changed and the writer will need an update, and set nFields = 0.
uint8_t direction
Definition: direction.H:48
A class for handling file names.
Definition: fileName.H:71
static const Enum< commsTypes > commsTypeNames
Enumerated names for the communication types.
Definition: UPstream.H:84
surfaceWriter()
Default construct.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition: label.H:164
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
bool empty() const
The surface to write is empty if the global number of faces is zero.
dictionary fieldLevel_
Field level to remove (on output)
void checkOpen() const
Verify that the outputPath_ has been set or FatalError.
virtual void open(const fileName &outputPath)
Open for output on specified path, using existing surface.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
#define defineSurfaceFieldMethods(ThisClass, Type)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
const meshedSurf & surface() const
Merge surfaces (if not upToDate) and return merged (parallel) or regular surface (non-parallel) ...
dictionary getFormatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Find "formatOptions" in a top-level dictionary. Extract and merge &#39;default&#39; + formatName values...
Definition: fileFormats.C:80
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1004
const meshedSurfRef & adjustSurface() const
Merge surfaces (if not upToDate) and return merged (parallel) or regular surface (non-parallel) and a...
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
point geometryCentre_
The centre of rotation (untranslate, translate)
runTimeSource setTime(sourceTimes[sourceTimeIndex], sourceTimeIndex)
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1184
static bool supportedType(const word &writeType)
True if New is likely to succeed for this writeType.
Definition: surfaceWriter.C:45
static bool canWriteType(const word &fileType, bool verbose=false)
Can this file format type be written via MeshedSurfaceProxy?
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition: UPstream.H:411
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
Definition: surfaceWriter.C:57
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.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
UPstream::commsTypes commType_
Communication type (for field merging)
Macros for easy insertion into run-time selection tables.
virtual void clear()
Reset origin and rotation to an identity coordinateSystem.
virtual ~surfaceWriter()
Destructor. Calls close()
Abstract definition of a meshed surface defined by faces and points.
Definition: meshedSurf.H:43
bool isPointData() const noexcept
Are the field data to be treated as point data?
A surfaceWriter that writes the geometry via the MeshedSurfaceProxy, but which does not support any f...
A Cartesian coordinate system.
Definition: cartesianCS.H:65
static scalar defaultMergeDim
The default merge dimension (1e-8)
virtual void clear()
Close any open output, remove association with a surface and expire the writer. The parallel flag rem...
virtual void setSurface(const meshedSurf &surf, bool parallel)
Change association with a surface, expire the writer with defined parallel/serial treatment...
Spatial transformation functions for primitive fields.
word timeName
Definition: getTimeIndex.H:3
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
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:62
virtual void endTime()
End a time-step.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const pointField & points
Implements a meshed surface by referencing another meshed surface or faces/points components...
Definition: meshedSurfRef.H:47
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual bool needsUpdate() const
Does the writer need an update (eg, lagging behind surface changes)
coordSystem::cartesian geometryTransform_
Local coordinate system transformation.
String literal.
Definition: keyType.H:82
scalar geometryScale_
Output geometry scaling after rotate/translate.
bool hasSurface() const
Writer is associated with a surface.
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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:770
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
tmp< Field< Type > > mergeFieldTemplate(const Field< Type > &fld) const
Gather (merge) fields with renumbering and shrinking for point data.
A proxy for writing MeshedSurface, UnsortedMeshedSurface and surfMesh to various file formats...
Definition: MeshedSurface.H:75
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary. ...
Definition: dictionary.C:533
virtual void close()
Finish output, performing any necessary cleanup.
tmp< Field< Type > > adjustFieldTemplate(const word &fieldName, const tmp< Field< Type >> &tfield) const
Apply refLevel and fieldScaling.
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))
dictionary fieldScale_
Field scaling (on output)
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name...
Definition: instant.H:53
A helper class for outputting values to Ostream.
Definition: ensightCells.H:43
virtual bool wroteData() const
Geometry or fields written since the last open?
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1037
bool readIfPresent(const word &key, const dictionary &dict, EnumType &val) const
Find an entry if present, and assign to T val.
Definition: EnumI.H:125
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool merge() const
Merge surfaces if they are not already upToDate (parallel) or simply mark the surface as being up-to-...
Reading is optional [identical to READ_IF_PRESENT].
Base class for surface writers.
bool verbose_
Additional output verbosity.
static void gather(const labelUList &offsets, const label comm, const ProcIDsContainer &procIDs, const UList< Type > &fld, List< Type > &allFld, const int tag=UPstream::msgType(), const UPstream::commsTypes=UPstream::commsTypes::nonBlocking)
Collect data in processor order on master (== procIDs[0]).
A class for managing temporary objects.
Definition: HashPtrTable.H:50
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
label & setComponent(label &val, const direction) noexcept
Non-const access to integer-type (has no components)
Definition: label.H:144
static wordHashSet writeTypes()
The file format types that can be written via MeshedSurfaceProxy.
Tensor of scalars, i.e. Tensor<scalar>.
Regular expression.
Definition: keyType.H:83
void unsetTime()
Clear the current time.
Inter-processor communications stream.
Definition: UPstream.H:62
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:80
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
virtual void beginTime(const Time &t)
Begin a time-step.
void clear()
Invalid by redirecting to null objects.
Namespace for OpenFOAM.
void setTime(const instant &inst)
Set the current time.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and a sub-dictionary) otherwise return nullptr...
Definition: dictionaryI.H:120
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133