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-2024 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::TryNew(const word& writeType)
88 {
89  // Without dictionary options
90  {
91  auto* ctorPtr = wordConstructorTable(writeType);
92 
93  if (ctorPtr)
94  {
95  return autoPtr<surfaceWriter>(ctorPtr());
96  }
97  }
98 
99  // Fallback to proxy writer...
100  return surfaceWriters::proxyWriter::TryNew(writeType);
101 }
102 
103 
106 (
107  const word& writeType,
108  const dictionary& writeOpts
109 )
110 {
111  // With dictionary options
112  {
113  auto* ctorPtr = wordDictConstructorTable(writeType);
114 
115  if (ctorPtr)
116  {
117  return autoPtr<surfaceWriter>(ctorPtr(writeOpts));
118  }
119  }
120 
121  // Without dictionary options
122  {
123  auto* ctorPtr = wordConstructorTable(writeType);
124 
125  if (ctorPtr)
126  {
127  return autoPtr<surfaceWriter>(ctorPtr());
128  }
129  }
131  // Fallback to proxy writer...
132  return surfaceWriters::proxyWriter::TryNew(writeType, writeOpts);
133 }
134 
135 
137 Foam::surfaceWriter::New(const word& writeType)
138 {
140  (
141  surfaceWriter::TryNew(writeType)
142  );
143 
144  if (!writer)
145  {
147  << "Unknown write type \"" << writeType << "\"\n\n"
148  << "Valid write types : "
149  << flatOutput(wordConstructorTablePtr_->sortedToc()) << nl
150  << "Valid proxy types : "
152  << exit(FatalError);
153  }
155  return writer;
156 }
157 
158 
161 (
162  const word& writeType,
163  const dictionary& writeOpts
164 )
165 {
167  (
168  surfaceWriter::TryNew(writeType, writeOpts)
169  );
170 
171  if (!writer)
172  {
174  << "Unknown write type \"" << writeType << "\"\n\n"
175  << "Valid write types : "
176  << flatOutput(wordConstructorTablePtr_->sortedToc()) << nl
177  << "Valid proxy types : "
179  << exit(FatalError);
180  }
182  return writer;
183 }
184 
185 
186 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
187 
189 :
190  surf_(),
191  mergedSurf_(),
192  adjustedSurf_(),
193  mergeDim_(defaultMergeDim),
194  geometryScale_(1),
195  geometryCentre_(Zero),
196  geometryTransform_(),
197  upToDate_(false),
198  wroteGeom_(false),
199  parallel_(true),
200  useTimeDir_(false),
201  isPointData_(false),
202  verbose_(false),
203  commType_(UPstream::commsTypes::scheduled),
204  nFields_(0),
205  currTime_(),
206  outputPath_(),
207  fieldLevel_(),
208  fieldScale_()
209 {
211 }
212 
213 
215 :
216  surfaceWriter()
217 {
218  options.readIfPresent("verbose", verbose_);
219 
220  UPstream::commsTypeNames.readIfPresent("commsType", options, commType_);
221 
222  geometryScale_ = 1;
225 
226  options.readIfPresent("scale", geometryScale_);
227 
228  // Optional cartesian coordinate system transform
229  const auto* dictptr = options.findDict("transform", keyType::LITERAL);
230 
231  if (dictptr)
232  {
233  dictptr->readIfPresent("rotationCentre", geometryCentre_);
234 
235  // 'origin' is optional within sub-dictionary
238  }
239 
240  fieldLevel_ = options.subOrEmptyDict("fieldLevel");
241  fieldScale_ = options.subOrEmptyDict("fieldScale");
242 
243  if (verbose_)
244  {
245  Info<< "Create surfaceWriter ("
246  << (this->isPointData() ? "point" : "face") << " data):"
247  << " commsType=" << UPstream::commsTypeNames[commType_] << endl;
248  }
249 }
250 
251 
252 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
253 
255 {
256  close();
257 }
258 
259 
260 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
262 void Foam::surfaceWriter::setTime(const instant& inst)
263 {
264  currTime_ = inst;
265 }
266 
268 void Foam::surfaceWriter::setTime(scalar timeValue)
269 {
270  currTime_ = instant(timeValue);
271 }
272 
273 
274 void Foam::surfaceWriter::setTime(scalar timeValue, const word& timeName)
275 {
276  currTime_.value() = timeValue;
277  currTime_.name() = timeName;
278 }
279 
280 
282 {
283  currTime_.value() = 0;
284  currTime_.name().clear();
285 }
286 
289 {
290  setTime(t.value(), t.timeName());
291 }
292 
294 void Foam::surfaceWriter::beginTime(const instant& inst)
295 {
296  setTime(inst);
297 }
298 
301 {
302  unsetTime();
303 }
304 
305 
306 void Foam::surfaceWriter::open(const fileName& outputPath)
307 {
308  outputPath_ = outputPath;
309  wroteGeom_ = false;
310 }
311 
312 
314 (
315  const meshedSurf& surf,
316  const fileName& outputPath,
317  bool parallel
318 )
319 {
320  close();
321  setSurface(surf, parallel);
322  open(outputPath);
323 }
324 
325 
327 (
328  const pointField& points,
329  const faceList& faces,
330  const fileName& outputPath,
331  bool parallel
332 )
333 {
334  close();
335  setSurface(points, faces, parallel);
336  open(outputPath);
337 }
338 
339 
341 (
342  const meshedSurf& surf,
343  const fileName& outputPath
344 )
345 {
346  close();
347  setSurface(surf, parallel_);
348  open(outputPath);
349 }
350 
351 
353 (
354  const pointField& points,
355  const faceList& faces,
356  const fileName& outputPath
357 )
358 {
359  close();
360  setSurface(points, faces, parallel_);
361  open(outputPath);
362 }
363 
364 
366 {
367  outputPath_.clear();
368  wroteGeom_ = false;
369 }
370 
371 
373 {
374  close();
375  expire();
376  surf_.clear();
377 }
378 
379 
381 (
382  const meshedSurf& surf,
383  bool parallel
384 )
385 {
386  expire();
387  surf_.reset(surf);
388  parallel_ = (parallel && UPstream::parRun());
389 }
390 
391 
393 (
394  const pointField& points,
395  const faceList& faces,
396  bool parallel
397 )
398 {
399  expire();
400  surf_.reset(points, faces);
401  parallel_ = (parallel && UPstream::parRun());
402 }
403 
404 
406 (
407  const meshedSurf& surf
408 )
409 {
410  setSurface(surf, parallel_);
411 }
412 
413 
415 (
416  const pointField& points,
417  const faceList& faces
418 )
419 {
420  setSurface(points, faces, parallel_);
421 }
422 
425 {
426  return !upToDate_;
427 }
428 
431 {
432  return wroteGeom_;
433 }
434 
435 
437 {
438  const bool changed = upToDate_;
439 
440  upToDate_ = false;
441  wroteGeom_ = false;
442  adjustedSurf_.clear();
443  mergedSurf_.clear();
444 
445  // Field count (nFields_) is a different type of accounting
446  // and is unaffected by geometry changes
447 
448  return changed;
449 }
450 
453 {
454  return surf_.good();
455 }
456 
457 
458 bool Foam::surfaceWriter::empty() const
459 {
460  const bool value = surf_.faces().empty();
461 
462  return (parallel_ ? returnReduceAnd(value) : value);
463 }
464 
465 
466 Foam::label Foam::surfaceWriter::size() const
467 {
468  const label value = surf_.faces().size();
469 
470  return (parallel_ ? returnReduce(value, sumOp<label>()) : value);
471 }
472 
473 
475 {
476  if (!is_open())
477  {
479  << type() << " : Attempted to write without a path" << nl
480  << exit(FatalError);
481  }
482 }
483 
484 
485 bool Foam::surfaceWriter::merge() const
486 {
487  bool changed = false;
488 
489  if (!upToDate_)
490  {
491  // Similar to expire
492  adjustedSurf_.clear();
493 
494  if (parallel_ && UPstream::parRun())
495  {
496  changed = mergedSurf_.merge(surf_, mergeDim_);
497  }
498  else
499  {
500  mergedSurf_.clear();
501  }
502 
503  if (changed)
504  {
505  wroteGeom_ = false;
506  }
507  }
508 
509  upToDate_ = true;
510  return changed;
511 }
512 
513 
515 {
516  merge();
517 
518  if (parallel_ && UPstream::parRun())
519  {
520  return mergedSurf_;
521  }
522 
523  return surf_;
524 }
525 
526 
528 {
529  if (!upToDate_)
530  {
531  adjustedSurf_.clear();
532  }
533 
534  if (!adjustedSurf_.good())
535  {
536  adjustedSurf_.reset(surface());
537 
538  tmp<pointField> tpts;
539 
540  if (geometryTransform_.good())
541  {
542  if (!geometryTransform_.R().is_identity())
543  {
544  if (magSqr(geometryCentre_) > ROOTVSMALL)
545  {
546  // Set centre of rotation,
547  // followed by forward transform (local -> global)
548  tpts =
549  geometryTransform_.globalPosition
550  (
551  adjustedSurf_.points0() - geometryCentre_
552  );
553 
554  // Unset centre of rotation
555  tpts.ref() += geometryCentre_;
556  }
557  else
558  {
559  // Forward transform (local -> global)
560  tpts =
561  geometryTransform_.globalPosition
562  (
563  adjustedSurf_.points0()
564  );
565  }
566  }
567  else if (magSqr(geometryTransform_.origin()) > ROOTVSMALL)
568  {
569  // Translate only (local -> global)
570  tpts = (adjustedSurf_.points0() + geometryTransform_.origin());
571  }
572  }
573 
574  adjustedSurf_.movePoints(tpts);
575  adjustedSurf_.scalePoints(geometryScale_);
576  }
577 
578  return adjustedSurf_;
579 }
580 
581 
582 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
583 
584 template<class Type>
586 (
587  const Field<Type>& fld
588 ) const
589 {
590  if (parallel_ && UPstream::parRun())
591  {
592  // Ensure geometry is also merged
593  merge();
594 
595  // Gather all values
596  auto tfield = tmp<Field<Type>>::New();
597  auto& allFld = tfield.ref();
598 
599  // Update any expired global index (as required)
600 
601  const globalIndex& globIndex =
602  (
603  this->isPointData()
604  ? mergedSurf_.pointGlobalIndex()
605  : mergedSurf_.faceGlobalIndex()
606  );
607 
608  globIndex.gather
609  (
610  fld,
611  allFld,
613  commType_,
615  );
616 
617  // Renumber (point data) to correspond to merged points
618  if
619  (
621  && this->isPointData()
622  && mergedSurf_.pointsMap().size()
623  )
624  {
625  inplaceReorder(mergedSurf_.pointsMap(), allFld);
626  allFld.resize(mergedSurf_.points().size());
627  }
628 
629  return tfield;
630  }
631 
632  // Mark that any geometry changes have been taken care of
633  upToDate_ = true;
635  return fld;
636 }
637 
638 
639 template<class Type>
641 (
642  const word& fieldName,
643  const tmp<Field<Type>>& tfield
644 ) const
645 {
646  if (verbose_)
647  {
648  Info<< "Writing field " << fieldName;
649  }
650 
651  tmp<Field<Type>> tadjusted;
652 
653  // Output scaling for the variable, but not for integer types
654  // which are typically ids etc.
655  if (!std::is_integral<Type>::value)
656  {
657  scalar value;
658 
659  // Remove *uniform* reference level
660  if
661  (
662  fieldLevel_.readIfPresent(fieldName, value, keyType::REGEX)
663  && !equal(value, 0)
664  )
665  {
666  // Could also detect brackets (...) and read accordingly
667  // or automatically scale by 1/sqrt(nComponents) instead ...
668 
669  Type refLevel;
670  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; ++cmpt)
671  {
672  setComponent(refLevel, cmpt) = value;
673  }
674 
675  if (verbose_)
676  {
677  Info<< " [level " << refLevel << ']';
678  }
679 
680  if (!tadjusted)
681  {
682  // Steal or clone
683  tadjusted.reset(tfield.ptr());
684  }
685 
686  // Remove offset level
687  tadjusted.ref() -= refLevel;
688  }
689 
690  // Apply scaling
691  if
692  (
693  fieldScale_.readIfPresent(fieldName, value, keyType::REGEX)
694  && !equal(value, 1)
695  )
696  {
697  if (verbose_)
698  {
699  Info<< " [scaling " << value << ']';
700  }
701 
702  if (!tadjusted)
703  {
704  // Steal or clone
705  tadjusted.reset(tfield.ptr());
706  }
707 
708  // Apply scaling
709  tadjusted.ref() *= value;
710  }
711 
712  // Rotate fields (vector and non-spherical tensors)
713  if
714  (
715  (is_vectorspace<Type>::value && pTraits<Type>::nComponents > 1)
716  && geometryTransform_.good()
717  && !geometryTransform_.R().is_identity()
718  )
719  {
720  if (!tadjusted)
721  {
722  // Steal or clone
723  tadjusted.reset(tfield.ptr());
724  }
725 
727  (
728  tadjusted.ref(),
729  geometryTransform_.R(),
730  tadjusted()
731  );
732  }
733  }
734 
735  return (tadjusted ? tadjusted : tfield);
736 }
737 
738 
739 #define defineSurfaceFieldMethods(ThisClass, Type) \
740  Foam::tmp<Foam::Field<Type>> \
741  ThisClass::mergeField(const Field<Type>& fld) const \
742  { \
743  return mergeFieldTemplate(fld); \
744  } \
745  \
746  Foam::tmp<Foam::Field<Type>> \
747  ThisClass::adjustField \
748  ( \
749  const word& fieldName, \
750  const tmp<Field<Type>>& tfield \
751  ) const \
752  { \
753  return adjustFieldTemplate(fieldName, tfield); \
754  }
755 
763 #undef defineSurfaceFieldMethod
764 
765 
766 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
767 
768 Foam::Ostream& Foam::operator<<
769 (
770  Ostream& os,
771  const InfoProxy<surfaceWriter>& iproxy
772 )
773 {
774  const auto& w = *iproxy;
775 
776  os << "surfaceWriter:"
777  << " upToDate: " << w.upToDate_
778  << " PointData: " << w.isPointData_
779  << " nFields: " << w.nFields_
780  << " time: " << w.currTime_
781  << " path: " << w.outputPath_ << endl;
782 
783  return os;
784 }
785 
786 
787 // ************************************************************************* //
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:46
A class for handling file names.
Definition: fileName.H:72
static const Enum< commsTypes > commsTypeNames
Enumerated names for the communication types.
Definition: UPstream.H:89
surfaceWriter()
Default construct.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
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:129
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:608
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)
bool readIfPresent(const word &key, const dictionary &dict, EnumType &val, const bool warnOnly=false) const
Find an entry if present, and assign to T val.
Definition: EnumI.H:111
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:50
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:531
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:1061
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:1252
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
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:421
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 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:61
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:56
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 a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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:521
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:44
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:1094
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
static autoPtr< surfaceWriter > TryNew(const word &writeType)
Optional select construct proxy writer.
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
Tensor of scalars, i.e. Tensor<scalar>.
static autoPtr< surfaceWriter > TryNew(const word &writeType)
Optional select construct surfaceWriter.
Definition: surfaceWriter.C:80
Regular expression.
Definition: keyType.H:83
void unsetTime()
Clear the current time.
Inter-processor communications stream.
Definition: UPstream.H:65
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
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 it is a dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127