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-------------------------------------------------------------------------------
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 "Time.H"
33#include "globalIndex.H"
34#include "coordinateRotation.H"
35#include "transformField.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
45}
46
48
49
50// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
51
53{
54 return
55 (
56 wordConstructorTablePtr_->found(writeType)
57 || wordDictConstructorTablePtr_->found(writeType)
59 );
60}
61
62
65{
66 // Constructors without dictionary options
67 auto* ctorPtr = wordConstructorTable(writeType);
68
69 if (!ctorPtr)
70 {
72 {
73 // Generally unknown, but handle via 'proxy' handler
75 (
76 new surfaceWriters::proxyWriter(writeType)
77 );
78 }
79
81 << "Unknown write type \"" << writeType << "\"\n\n"
82 << "Valid write types : "
83 << flatOutput(wordConstructorTablePtr_->sortedToc()) << nl
84 << "Valid proxy types : "
86 << exit(FatalError);
87 }
88
89 return autoPtr<surfaceWriter>(ctorPtr());
90}
91
92
95(
96 const word& writeType,
97 const dictionary& writeOpts
98)
99{
100 // Constructors with dictionary options
101 {
102 auto* ctorPtr = wordDictConstructorTable(writeType);
103
104 if (ctorPtr)
105 {
106 return autoPtr<surfaceWriter>(ctorPtr(writeOpts));
107 }
108 }
109
110
111 // Constructors without dictionary options
112 auto* ctorPtr = wordConstructorTable(writeType);
113
114 if (!ctorPtr)
115 {
117 {
118 // Generally unknown, but handle via 'proxy' handler
120 (
121 new surfaceWriters::proxyWriter(writeType, writeOpts)
122 );
123 }
124
126 << "Unknown write type \"" << writeType << "\"\n\n"
127 << "Valid write types : "
128 << wordConstructorTablePtr_->sortedToc() << nl
129 << "Valid proxy types : "
131 << exit(FatalError);
132 }
133
134 return autoPtr<surfaceWriter>(ctorPtr());
135}
136
137
138// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
139
141:
142 surf_(),
143 mergedSurf_(),
144 adjustedSurf_(),
145 mergeDim_(defaultMergeDim),
146 geometryScale_(1),
147 geometryTransform_(),
148 upToDate_(false),
149 wroteGeom_(false),
150 parallel_(true),
151 useTimeDir_(false),
152 isPointData_(false),
153 verbose_(false),
154 nFields_(0),
155 currTime_(),
156 outputPath_(),
157 fieldLevel_(),
158 fieldScale_()
159{
161}
162
163
165:
167{
168 options.readIfPresent("verbose", verbose_);
169
170 geometryScale_ = 1;
172
173 options.readIfPresent("scale", geometryScale_);
174
175 const dictionary* dictptr;
176
177 // Optional cartesian coordinate system transform
178 if ((dictptr = options.findDict("transform", keyType::LITERAL))!= nullptr)
179 {
181 }
182
183 fieldLevel_ = options.subOrEmptyDict("fieldLevel");
184 fieldScale_ = options.subOrEmptyDict("fieldScale");
185}
186
187
188// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
189
191{
192 close();
193}
194
195
196// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
197
199{
200 currTime_ = inst;
201}
202
203
204void Foam::surfaceWriter::setTime(scalar timeValue)
205{
206 currTime_ = instant(timeValue);
207}
208
209
210void Foam::surfaceWriter::setTime(scalar timeValue, const word& timeName)
211{
212 currTime_.value() = timeValue;
213 currTime_.name() = timeName;
214}
215
216
218{
219 currTime_.value() = 0;
220 currTime_.name().clear();
221}
222
223
225{
226 setTime(t.value(), t.timeName());
227}
228
229
231{
232 setTime(inst);
233}
234
235
237{
238 unsetTime();
239}
240
241
242void Foam::surfaceWriter::open(const fileName& outputPath)
243{
244 outputPath_ = outputPath;
245 wroteGeom_ = false;
246}
247
248
250(
251 const meshedSurf& surf,
252 const fileName& outputPath,
253 bool parallel
254)
255{
256 close();
257 setSurface(surf, parallel);
258 open(outputPath);
259}
260
261
263(
264 const pointField& points,
265 const faceList& faces,
266 const fileName& outputPath,
267 bool parallel
268)
269{
270 close();
271 setSurface(points, faces, parallel);
272 open(outputPath);
273}
274
275
277(
278 const meshedSurf& surf,
279 const fileName& outputPath
280)
281{
282 close();
283 setSurface(surf, parallel_);
284 open(outputPath);
285}
286
287
289(
290 const pointField& points,
291 const faceList& faces,
292 const fileName& outputPath
293)
294{
295 close();
296 setSurface(points, faces, parallel_);
297 open(outputPath);
298}
299
300
302{
303 outputPath_.clear();
304 wroteGeom_ = false;
305}
306
307
309{
310 close();
311 expire();
312 surf_.clear();
313}
314
315
317(
318 const meshedSurf& surf,
319 bool parallel
320)
321{
322 expire();
323 surf_.reset(surf);
324 parallel_ = (parallel && Pstream::parRun());
325}
326
327
329(
330 const pointField& points,
331 const faceList& faces,
332 bool parallel
333)
334{
335 expire();
336 surf_.reset(points, faces);
337 parallel_ = (parallel && Pstream::parRun());
338}
339
340
342(
343 const meshedSurf& surf
344)
345{
346 setSurface(surf, parallel_);
347}
348
349
351(
352 const pointField& points,
353 const faceList& faces
354)
355{
356 setSurface(points, faces, parallel_);
357}
358
359
361{
362 return !upToDate_;
363}
364
365
367{
368 return wroteGeom_;
369}
370
371
373{
374 const bool changed = upToDate_;
375
376 upToDate_ = false;
377 wroteGeom_ = false;
378 adjustedSurf_.clear();
379 mergedSurf_.clear();
380
381 // Field count (nFields_) is a different type of accounting
382 // and is unaffected by geometry changes
383
384 return changed;
385}
386
387
389{
390 return surf_.valid();
391}
392
393
395{
396 const bool value = surf_.faces().empty();
397
398 return (parallel_ ? returnReduce(value, andOp<bool>()) : value);
399}
400
401
402Foam::label Foam::surfaceWriter::size() const
403{
404 const label value = surf_.faces().size();
405
406 return (parallel_ ? returnReduce(value, sumOp<label>()) : value);
407}
408
409
411{
412 if (!is_open())
413 {
415 << type() << " : Attempted to write without a path" << nl
416 << exit(FatalError);
417 }
418}
419
420
422{
423 bool changed = false;
424
425 if (!upToDate_)
426 {
427 adjustedSurf_.clear();
428
429 if (parallel_ && Pstream::parRun())
430 {
431 changed = mergedSurf_.merge(surf_, mergeDim_);
432 }
433 else
434 {
435 mergedSurf_.clear();
436 }
437 }
438
439 upToDate_ = true;
440
441 if (changed)
442 {
443 wroteGeom_ = false;
444 }
445
446 return changed;
447}
448
449
451{
452 merge();
453
454 if (parallel_ && Pstream::parRun())
455 {
456 return mergedSurf_;
457 }
458
459 return surf_;
460}
461
462
464{
465 if (!upToDate_)
466 {
467 adjustedSurf_.clear();
468 }
469
470 if (!adjustedSurf_.valid())
471 {
472 adjustedSurf_.reset(surface());
473
474 if
475 (
476 geometryTransform_.valid()
477 &&
478 (
479 (magSqr(geometryTransform_.origin()) > ROOTVSMALL)
480 || !geometryTransform_.R().is_identity()
481 )
482 )
483 {
484 // Forward transform
485 adjustedSurf_.movePoints
486 (
487 geometryTransform_.globalPosition(adjustedSurf_.points0())
488 );
489 }
490
491 adjustedSurf_.scalePoints(geometryScale_);
492 }
493
494 return adjustedSurf_;
495}
496
497
498// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
499
500template<class Type>
502(
503 const Field<Type>& fld
504) const
505{
506 if (parallel_ && Pstream::parRun())
507 {
508 // Ensure geometry is also merged
509 merge();
510
511 // Gather all values
512 auto tfield = tmp<Field<Type>>::New();
513 auto& allFld = tfield.ref();
514
515 globalIndex::gatherOp(fld, allFld);
516
517 // Renumber (point data) to correspond to merged points
518 if
519 (
521 && this->isPointData()
522 && mergedSurf_.pointsMap().size()
523 )
524 {
525 inplaceReorder(mergedSurf_.pointsMap(), allFld);
526 allFld.resize(mergedSurf_.points().size());
527 }
528
529 return tfield;
530 }
531
532 // Mark that any geometry changes have been taken care of
533 upToDate_ = true;
534
535 return fld;
536}
537
538
539template<class Type>
541(
542 const word& fieldName,
543 const tmp<Field<Type>>& tfield
544) const
545{
546 if (verbose_)
547 {
548 Info<< "Writing field " << fieldName;
549 }
550
551 tmp<Field<Type>> tadjusted;
552
553 // Output scaling for the variable, but not for integer types
554 // which are typically ids etc.
555 if (!std::is_integral<Type>::value)
556 {
557 scalar value;
558
559 // Remove *uniform* reference level
560 if
561 (
562 fieldLevel_.readIfPresent(fieldName, value, keyType::REGEX)
563 && !equal(value, 0)
564 )
565 {
566 // Could also detect brackets (...) and read accordingly
567 // or automatically scale by 1/sqrt(nComponents) instead ...
568
569 Type refLevel;
570 for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; ++cmpt)
571 {
572 setComponent(refLevel, cmpt) = value;
573 }
574
575 if (verbose_)
576 {
577 Info<< " [level " << refLevel << ']';
578 }
579
580 if (!tadjusted)
581 {
582 // Steal or clone
583 tadjusted.reset(tfield.ptr());
584 }
585
586 // Remove offset level
587 tadjusted.ref() -= refLevel;
588 }
589
590 // Apply scaling
591 if
592 (
593 fieldScale_.readIfPresent(fieldName, value, keyType::REGEX)
594 && !equal(value, 1)
595 )
596 {
597 if (verbose_)
598 {
599 Info<< " [scaling " << value << ']';
600 }
601
602 if (!tadjusted)
603 {
604 // Steal or clone
605 tadjusted.reset(tfield.ptr());
606 }
607
608 // Apply scaling
609 tadjusted.ref() *= value;
610 }
611
612 // Rotate fields (vector and non-spherical tensors)
613 if
614 (
616 && geometryTransform_.valid()
617 && !geometryTransform_.R().is_identity()
618 )
619 {
620 if (!tadjusted)
621 {
622 // Steal or clone
623 tadjusted.reset(tfield.ptr());
624 }
625
627 (
628 tadjusted.ref(),
629 geometryTransform_.R(),
630 tadjusted()
631 );
632 }
633 }
634
635 return (tadjusted ? tadjusted : tfield);
636}
637
638
639#define defineSurfaceFieldMethods(ThisClass, Type) \
640 Foam::tmp<Foam::Field<Type>> \
641 ThisClass::mergeField(const Field<Type>& fld) const \
642 { \
643 return mergeFieldTemplate(fld); \
644 } \
645 \
646 Foam::tmp<Foam::Field<Type>> \
647 ThisClass::adjustField \
648 ( \
649 const word& fieldName, \
650 const tmp<Field<Type>>& tfield \
651 ) const \
652 { \
653 return adjustFieldTemplate(fieldName, tfield); \
654 }
655
662
663#undef defineSurfaceFieldMethod
664
665
666// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
667
668Foam::Ostream& Foam::operator<<
669(
670 Ostream& os,
672)
673{
674 const surfaceWriter& w = ip.t_;
675
676 os << "surfaceWriter:"
677 << " upToDate: " << w.upToDate_
678 << " PointData: " << w.isPointData_
679 << " nFields: " << w.nFields_
680 << " time: " << w.currTime_
681 << " path: " << w.outputPath_ << endl;
682
683 return os;
684}
685
686
687// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:52
A proxy for writing MeshedSurface, UnsortedMeshedSurface and surfMesh to various file formats.
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:62
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A Cartesian coordinate system.
Definition: cartesianCS.H:72
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:126
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
const Type & value() const
Return const reference to value.
A class for handling file names.
Definition: fileName.H:76
static void gatherOp(const UList< Type > &sendData, List< Type > &allData, const int tag=UPstream::msgType(), const UPstream::commsTypes=UPstream::commsTypes::nonBlocking, const label comm=UPstream::worldComm)
Collect data in processor order on master (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:81
@ REGEX
Regular expression.
Definition: keyType.H:82
Implements a meshed surface by referencing another meshed surface or faces/points components.
Definition: meshedSurfRef.H:56
void clear()
Invalid by redirecting to null objects.
Abstract definition of a meshed surface defined by faces and points.
Definition: meshedSurf.H:50
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
splitCell * master() const
Definition: splitCell.H:113
Base class for surface writers.
surfaceWriter()
Default construct.
coordSystem::cartesian geometryTransform_
Local coordinate system transformation.
virtual void endTime()
End a time-step.
scalar geometryScale_
Output geometry scaling after rotate/translate.
bool isPointData_
Is point vs cell data.
const meshedSurf & surface() const
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)
virtual void close()
Finish output, performing any necessary cleanup.
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?
dictionary fieldLevel_
Field level to remove (on output)
static bool supportedType(const word &writeType)
True if New is likely to succeed for this writeType.
Definition: surfaceWriter.C:52
void setTime(const instant &inst)
Set the current time.
virtual bool merge() const
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()
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
virtual void clear()
virtual bool wroteData() const
Geometry or fields written since the last open?
fileName outputPath_
The full output directory and file (surface) name.
A surfaceWriter that writes the geometry via the MeshedSurfaceProxy, but which does not support any f...
A class for managing temporary objects.
Definition: tmp.H:65
void reset(tmp< T > &&other) noexcept
Clear existing and transfer ownership.
Definition: tmpI.H:314
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
runTimeSource setTime(sourceTimes[sourceTimeIndex], sourceTimeIndex)
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
label & setComponent(label &val, const direction) noexcept
Non-const access to integer-type (has no components)
Definition: label.H:123
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
bool equal(const T &s1, const T &s2)
Compare two values for equality.
Definition: doubleFloat.H:46
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
uint8_t direction
Definition: direction.H:56
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
volScalarField & e
Definition: createFields.H:11
static constexpr char open
Definition: FlatOutput.H:73
#define defineSurfaceFieldMethods(ThisClass, Type)
Spatial transformation functions for primitive fields.