Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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"
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
46
48
49
50// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
51
52bool 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
81{
82 return fileFormats::getFormatOptions(dict, surfDict, formatName, entryName);
83}
84
85
87Foam::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 }
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
137Foam::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);
154
155 return writer;
156}
157
158
161(
162 const word& writeType,
163 const dictionary& writeOpts
164)
165{
166 autoPtr<surfaceWriter> writer
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 gatherv_(false),
205 nFields_(0),
206 currTime_(),
207 outputPath_(),
210{
212}
213
214
216:
218{
219 options.readIfPresent("verbose", verbose_);
220
221 UPstream::commsTypeNames.readIfPresent("commsType", options, commType_);
222 gatherv_ = false;
223 options.readIfPresent("gatherv", gatherv_);
224
225 geometryScale_ = 1;
228
229 options.readIfPresent("scale", geometryScale_);
230
231 // Optional cartesian coordinate system transform
232 const auto* dictptr = options.findDict("transform", keyType::LITERAL);
233
234 if (dictptr)
235 {
236 dictptr->readIfPresent("rotationCentre", geometryCentre_);
237
238 // 'origin' is optional within sub-dictionary
241 }
242
243 fieldLevel_ = options.subOrEmptyDict("fieldLevel");
244 fieldScale_ = options.subOrEmptyDict("fieldScale");
245
246 if (verbose_)
247 {
248 Info<< "Create surfaceWriter ("
249 << (this->isPointData() ? "point" : "face") << " data):"
250 << " commsType=";
251
252 if (UPstream::parRun())
253 {
254 if (gatherv_) Info<< "gatherv+";
256 }
257 else
258 {
259 Info<< "serial";
260 }
261
263 }
264}
265
266
267// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
268
271 close();
272}
273
274
275// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
278{
279 currTime_ = inst;
280}
281
287
288
290{
291 currTime_.value() = timeValue;
292 currTime_.name() = timeName;
293}
294
295
297{
298 currTime_.value() = 0;
299 currTime_.name().clear();
300}
301
304{
305 setTime(t.value(), t.timeName());
306}
307
310{
311 setTime(inst);
312}
313
316{
317 unsetTime();
318}
319
320
321void Foam::surfaceWriter::open(const fileName& outputPath)
322{
323 outputPath_ = outputPath;
324 wroteGeom_ = false;
325}
326
327
329(
330 const meshedSurf& surf,
331 const fileName& outputPath,
332 bool parallel
333)
335 close();
336 setSurface(surf, parallel);
337 open(outputPath);
338}
339
340
342(
343 const pointField& points,
344 const faceList& faces,
345 const fileName& outputPath,
346 bool parallel
347)
349 close();
350 setSurface(points, faces, parallel);
351 open(outputPath);
352}
353
354
356(
357 const meshedSurf& surf,
358 const fileName& outputPath
359)
361 close();
362 setSurface(surf, parallel_);
363 open(outputPath);
364}
365
366
368(
369 const pointField& points,
370 const faceList& faces,
371 const fileName& outputPath
372)
374 close();
375 setSurface(points, faces, parallel_);
376 open(outputPath);
377}
378
379
381{
382 outputPath_.clear();
383 wroteGeom_ = false;
384}
385
386
389 close();
390 expire();
391 surf_.clear();
392}
393
394
396(
397 const meshedSurf& surf,
398 bool parallel
399)
401 expire();
402 surf_.reset(surf);
403 parallel_ = (parallel && UPstream::parRun());
404}
405
406
408(
409 const pointField& points,
410 const faceList& faces,
411 bool parallel
412)
414 expire();
415 surf_.reset(points, faces);
416 parallel_ = (parallel && UPstream::parRun());
417}
418
419
421(
422 const meshedSurf& surf
423)
424{
425 setSurface(surf, parallel_);
426}
427
428
430(
431 const pointField& points,
432 const faceList& faces
433)
434{
435 setSurface(points, faces, parallel_);
436}
437
440{
441 return !upToDate_;
442}
443
446{
447 return wroteGeom_;
448}
449
450
452{
453 const bool changed = upToDate_;
454
455 upToDate_ = false;
456 wroteGeom_ = false;
457 adjustedSurf_.clear();
458 mergedSurf_.clear();
459
460 // Field count (nFields_) is a different type of accounting
461 // and is unaffected by geometry changes
462
463 return changed;
464}
465
468{
469 return surf_.good();
470}
471
472
475 const bool value = surf_.faces().empty();
476
477 return (parallel_ ? returnReduceAnd(value) : value);
478}
479
480
481Foam::label Foam::surfaceWriter::size() const
483 const label value = surf_.faces().size();
484
485 return (parallel_ ? returnReduce(value, sumOp<label>()) : value);
486}
487
488
490{
491 if (!is_open())
492 {
494 << type() << " : Attempted to write without a path" << nl
495 << exit(FatalError);
496 }
497}
498
499
501{
502 bool changed = false;
503
504 if (!upToDate_)
505 {
506 // Similar to expire
507 adjustedSurf_.clear();
508
509 if (parallel_ && UPstream::parRun())
510 {
511 changed = mergedSurf_.merge(surf_, mergeDim_);
512 }
513 else
514 {
515 mergedSurf_.clear();
516 }
517
518 if (changed)
519 {
520 wroteGeom_ = false;
521 }
523
524 upToDate_ = true;
525 return changed;
526}
527
528
530{
531 merge();
532
533 if (parallel_ && UPstream::parRun())
534 {
536 }
537
538 return surf_;
539}
540
541
543{
544 if (!upToDate_)
545 {
546 adjustedSurf_.clear();
547 }
548
549 if (!adjustedSurf_.good())
550 {
551 adjustedSurf_.reset(surface());
552
553 tmp<pointField> tpts;
554
555 if (geometryTransform_.good())
556 {
557 if (!geometryTransform_.R().is_identity())
558 {
559 if (magSqr(geometryCentre_) > ROOTVSMALL)
560 {
561 // Set centre of rotation,
562 // followed by forward transform (local -> global)
563 tpts =
564 geometryTransform_.globalPosition
565 (
566 adjustedSurf_.points0() - geometryCentre_
567 );
568
569 // Unset centre of rotation
570 tpts.ref() += geometryCentre_;
571 }
572 else
573 {
574 // Forward transform (local -> global)
575 tpts =
576 geometryTransform_.globalPosition
577 (
578 adjustedSurf_.points0()
579 );
580 }
581 }
582 else if (magSqr(geometryTransform_.origin()) > ROOTVSMALL)
583 {
584 // Translate only (local -> global)
585 tpts = (adjustedSurf_.points0() + geometryTransform_.origin());
586 }
587 }
588
589 adjustedSurf_.movePoints(tpts);
590 adjustedSurf_.scalePoints(geometryScale_);
591 }
592
594}
595
596
597// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
598
599template<class Type>
601(
602 const Field<Type>& fld
603) const
604{
605 if (parallel_ && UPstream::parRun())
606 {
607 // Ensure geometry is also merged
608 merge();
609
610 // Gather all values
611 auto tfield = tmp<Field<Type>>::New();
612 auto& allFld = tfield.ref();
613
614 // Update any expired global index (as required)
615
616 const globalIndex& globIndex =
617 (
618 this->isPointData()
619 ? mergedSurf_.pointGlobalIndex()
620 : mergedSurf_.faceGlobalIndex()
621 );
622
623 if (gatherv_)
624 {
625 globIndex.mpiGather
626 (
627 fld,
628 allFld,
630 // For fallback:
631 commType_,
633 );
634 }
635 else
636 {
637 globIndex.gather
638 (
639 fld,
640 allFld,
642 commType_,
644 );
645 }
646
647 // Renumber (point data) to correspond to merged points
648 if
649 (
651 && this->isPointData()
652 && mergedSurf_.pointsMap().size()
653 )
654 {
655 inplaceReorder(mergedSurf_.pointsMap(), allFld);
656 allFld.resize(mergedSurf_.points().size());
657 }
658
659 // Extended debugging. Limit to master:
660 #if 0
661 if (UPstream::master())
662 {
663 Info<< "merged List<" << pTraits<Type>::typeName << "> : ";
664 allFld.writeList(Info) << endl;
665 }
666 #endif
667
668 return tfield;
669 }
670
671 // Mark that any geometry changes have been taken care of
672 upToDate_ = true;
673
674 return fld;
675}
676
677
678template<class Type>
680(
681 const word& fieldName,
682 const tmp<Field<Type>>& tfield
683) const
684{
685 if (verbose_)
686 {
687 Info<< "Writing field " << fieldName;
688 }
689
690 tmp<Field<Type>> tadjusted;
691
692 // Output scaling for the variable, but not for integer types
693 // which are typically ids etc.
694 if constexpr (std::is_integral_v<Type>)
695 {
696 return tfield;
697 }
698 else
699 {
700 scalar value;
701
702 // Remove *uniform* reference level
703 if
704 (
705 fieldLevel_.readIfPresent(fieldName, value, keyType::REGEX)
706 && !equal(value, 0)
707 )
708 {
709 // Could also detect brackets (...) and read accordingly
710 // or automatically scale by 1/sqrt(nComponents) instead ...
711
712 Type refLevel;
713 if constexpr (is_vectorspace_v<Type>)
714 {
715 refLevel.fill(value);
716 }
717 else
718 {
719 refLevel = value;
720 }
721
722 if (verbose_)
723 {
724 Info<< " [level " << refLevel << ']';
725 }
726
727 if (!tadjusted)
728 {
729 // Steal or clone
730 tadjusted.reset(tfield.ptr());
731 }
732
733 // Remove offset level
734 tadjusted.ref() -= refLevel;
735 }
736
737 // Apply scaling
738 if
739 (
740 fieldScale_.readIfPresent(fieldName, value, keyType::REGEX)
741 && !equal(value, 1)
742 )
743 {
744 if (verbose_)
745 {
746 Info<< " [scaling " << value << ']';
747 }
748
749 if (!tadjusted)
750 {
751 // Steal or clone
752 tadjusted.reset(tfield.ptr());
753 }
754
755 // Apply scaling
756 tadjusted.ref() *= value;
757 }
758 }
759
760
761 // Rotate fields (vector and non-spherical tensors)
763 {
764 if (geometryTransform_.good() && !geometryTransform_.R().is_identity())
765 {
766 if (!tadjusted)
767 {
768 // Steal or clone
769 tadjusted.reset(tfield.ptr());
770 }
771
773 (
774 tadjusted.ref(),
775 geometryTransform_.R(),
776 tadjusted()
777 );
779 }
780
781 return (tadjusted ? tadjusted : tfield);
782}
783
784
785#define defineSurfaceFieldMethods(ThisClass, Type) \
786 Foam::tmp<Foam::Field<Type>> \
787 ThisClass::mergeField(const Field<Type>& fld) const \
788 { \
789 return mergeFieldTemplate(fld); \
790 } \
791 \
792 Foam::tmp<Foam::Field<Type>> \
793 ThisClass::adjustField \
794 ( \
795 const word& fieldName, \
796 const tmp<Field<Type>>& tfield \
797 ) const \
798 { \
799 return adjustFieldTemplate(fieldName, tfield); \
801
808
809#undef defineSurfaceFieldMethod
810
811
812// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
813
814Foam::Ostream& Foam::operator<<
815(
816 Ostream& os,
817 const InfoProxy<surfaceWriter>& iproxy
818)
819{
820 const auto& w = *iproxy;
821
822 os << "surfaceWriter:"
823 << " upToDate: " << w.upToDate_
824 << " PointData: " << w.isPointData_
825 << " nFields: " << w.nFields_
826 << " time: " << w.currTime_
827 << " path: " << w.outputPath_ << endl;
828
829 return os;
830}
831
832
833// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
@ LAZY_READ
Reading is optional [identical to READ_IF_PRESENT].
A helper class for outputting values to Ostream.
Definition InfoProxy.H:49
static bool canWriteType(const word &fileType, bool verbose=false)
Can this file format type be written via MeshedSurfaceProxy?
static wordHashSet writeTypes()
The file format types that can be written via MeshedSurfaceProxy.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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
Inter-processor communications stream.
Definition UPstream.H:69
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static const Enum< commsTypes > commsTypeNames
Enumerated names for the communication types.
Definition UPstream.H:92
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition UPstream.H:1069
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A Cartesian coordinate system.
Definition cartesianCS.H:68
virtual void clear()
Reset origin and rotation to an identity coordinateSystem.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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
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...
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...
const Type & value() const noexcept
Return const reference to value.
A class for handling file names.
Definition fileName.H:75
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
static void gather(const labelUList &offsets, const label comm, const ProcIDsContainer &procIDs, const UList< Type > &fld, UList< Type > &allFld, const int tag=UPstream::msgType(), UPstream::commsTypes commsType=UPstream::commsTypes::nonBlocking)
Collect data in processor order on master (== procIDs[0]).
void mpiGather(const UList< Type > &sendData, OutputContainer &allData, const label comm=UPstream::worldComm, UPstream::commsTypes commsType=UPstream::commsTypes::nonBlocking, const int tag=UPstream::msgType()) const
Use MPI_Gatherv call for contiguous data when possible (in serial: performs a simple copy).
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name.
Definition instant.H:56
@ LITERAL
String literal.
Definition keyType.H:82
@ REGEX
Regular expression.
Definition keyType.H:83
Implements a meshed surface by referencing another meshed surface or faces/points components.
Abstract definition of a meshed surface defined by faces and points.
Definition meshedSurf.H:44
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Base class for surface writers.
surfaceWriter()
Default construct.
coordSystem::cartesian geometryTransform_
Local coordinate system transformation.
virtual void open(const fileName &outputPath)
Open for output on specified path, using existing surface.
meshedSurfRef surf_
Reference to surface or surface components.
virtual void endTime()
End a time-step.
scalar geometryScale_
Output geometry scaling after rotate/translate.
bool useTimeDir_
Insert additional time sub-directory in the output path.
bool gatherv_
Prefer MPI gatherv intrinsic (for field merging) [experimental].
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
bool wroteGeom_
Track if geometry has been written since the last open.
bool isPointData_
Is point vs cell data.
meshedSurfRef adjustedSurf_
The surface after point coordinate transforms and scaling.
const meshedSurf & surface() const
Merge surfaces (if not upToDate) and return merged (parallel) or regular surface (non-parallel).
bool isPointData() const noexcept
Are the field data to be treated as point data?
const word & timeName() const
The current time value/name.
void checkOpen() const
Verify that the outputPath_ has been set or FatalError.
void unsetTime()
Clear the current time.
virtual void beginTime(const Time &t)
Begin a time-step.
label size() const
The global number of faces for the associated surface.
label nFields_
The number of fields.
static scalar defaultMergeDim
The default merge dimension (1e-8).
virtual ~surfaceWriter()
Destructor. Calls close().
virtual void setSurface(const meshedSurf &surf, bool parallel)
Change association with a surface, expire the writer with defined parallel/serial treatment.
static autoPtr< surfaceWriter > TryNew(const word &writeType)
Optional select construct surfaceWriter.
virtual void close()
Finish output, performing any necessary cleanup.
bool parallel_
Writing in parallel (via master).
bool empty() const
The surface to write is empty if the global number of faces is zero.
bool upToDate_
The topology/surface is up-to-date?
UPstream::commsTypes commType_
Communication type (for field merging).
dictionary fieldLevel_
Field level to remove (on output).
static bool supportedType(const word &writeType)
True if New is likely to succeed for this writeType.
bool is_open() const noexcept
Test if outputPath has been set.
void setTime(const instant &inst)
Set the current time.
virtual bool merge() const
Merge surfaces if they are not already upToDate (parallel) or simply mark the surface as being up-to-...
mergedSurf mergedSurf_
Surface after merging (parallel).
bool verbose_
Additional output verbosity.
tmp< Field< Type > > adjustFieldTemplate(const word &fieldName, const tmp< Field< Type > > &tfield) const
Apply refLevel and fieldScaling.
dictionary fieldScale_
Field scaling (on output).
virtual bool expire()
Mark that surface changed and the writer will need an update, and set nFields = 0.
virtual bool needsUpdate() const
Does the writer need an update (eg, lagging behind surface changes).
tmp< Field< Type > > mergeFieldTemplate(const Field< Type > &fld) const
Gather (merge) fields with renumbering and shrinking for point data.
instant currTime_
The current time value/name.
bool hasSurface() const
Writer is associated with a surface.
const meshedSurfRef & adjustSurface() const
Merge surfaces (if not upToDate) and return merged (parallel) or regular surface (non-parallel) and a...
virtual void clear()
Close any open output, remove association with a surface and expire the writer. The parallel flag rem...
point geometryCentre_
The centre of rotation (untranslate, translate).
scalar timeValue() const
The current time value/name.
virtual bool wroteData() const
Geometry or fields written since the last open?
fileName outputPath_
The full output directory and file (surface) name.
scalar mergeDim_
Dimension for merging.
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
static autoPtr< surfaceWriter > TryNew(const word &writeType)
Optional select construct proxy writer.
A class for managing temporary objects.
Definition tmp.H:75
void reset(tmp< T > &&other) noexcept
Clear existing and transfer ownership.
Definition tmpI.H:338
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
runTimeSource setTime(sourceTimes[sourceTimeIndex], sourceTimeIndex)
const wordList surface
Standard surface field types (scalar, vector, tensor, etc).
dictionary getFormatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Find "formatOptions" in a top-level dictionary. Extract and merge 'default' + formatName values.
Definition fileFormats.C:80
Namespace for OpenFOAM.
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.
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition label.H:180
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Tensor< scalar > tensor
Definition symmTensor.H:57
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr bool is_vectorspace_v
The is_vectorspace value of Type.
Definition pTraits.H:244
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
constexpr bool is_rotational_vectorspace_v
The is_rotational_vectorspace value of Type.
Definition pTraits.H:251
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
SphericalTensor< scalar > sphericalTensor
SphericalTensor of scalars, i.e. SphericalTensor<scalar>.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & e
#define defineSurfaceFieldMethods(ThisClass, Type)
Spatial transformation functions for primitive fields.