gltfCoordSetWriter.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) 2021-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 "gltfCoordSetWriter.H"
29#include "coordSet.H"
30#include "fileName.H"
31#include "OFstream.H"
32#include "OSspecific.H"
33#include "foamGltfScene.H"
34#include "foamGltfSceneWriter.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
42namespace coordSetWriters
43{
47}
48}
49
52({
53 // No naming for NONE
54 { fieldOption::UNIFORM, "uniform" },
55 { fieldOption::FIELD, "field" },
56});
57
58
59// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
60
61namespace Foam
62{
63
64template<class Type>
66(
67 const colourTable& colours,
68 const Field<Type>& field,
69 const scalar boundMin,
70 const scalar boundMax
71)
72{
73 const label boundDelta = (boundMax - boundMin + ROOTVSMALL);
74
75 auto tresult = tmp<vectorField>::New(field.size());
76 auto& result = tresult.ref();
77
78 forAll(field, i)
79 {
80 const Type& val = field[i];
81
82 const scalar f =
83 (
85 ? scalar(component(val, 0))
86 : scalar(mag(val))
87 );
88
89 // 0-1 clipped by value()
90 result[i] = colours.value((f - boundMin)/boundDelta);
91 }
92
93 return tresult;
94}
95
96
97template<class Type>
99(
100 const dictionary& dict,
101 const colourTable& colours,
102 const Field<Type>& field
103)
104{
105 scalar refValue(0);
106 scalarMinMax valLimits;
107
109 {
110 MinMax<Type> scanned(minMax(field));
111
112 refValue = scalar(component(field[0], 0));
113 valLimits.min() = scalar(component(scanned.min(), 0));
114 valLimits.max() = scalar(component(scanned.max(), 0));
115 }
116 else
117 {
118 // Use mag() for multiple components
119 refValue = mag(field[0]);
120 valLimits = minMaxMag(field);
121 }
122
123 dict.readIfPresent("min", valLimits.min());
124 dict.readIfPresent("max", valLimits.max());
125
126 const scalar fraction =
127 (
128 (refValue - valLimits.min())
129 / (valLimits.max() - valLimits.min() + ROOTVSMALL)
130 );
131
132 // 0-1 clipped by value()
133 return colours.value(fraction);
134}
135
136
137} // End namespace Foam
138
139
140// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
141
142Foam::word Foam::coordSetWriters::gltfWriter::getColourMap
143(
144 const dictionary& dict
145) const
146{
148 dict.readIfPresent("colourMap", colourMap);
149
150 return colourMap;
151}
152
153
154const Foam::colourTable& Foam::coordSetWriters::gltfWriter::getColourTable
155(
156 const dictionary& dict
157) const
158{
159 return colourTable::ref(getColourMap(dict));
160}
161
162
163Foam::scalarMinMax Foam::coordSetWriters::gltfWriter::getFieldLimits
164(
165 const word& fieldName
166) const
167{
168 const dictionary fieldDict = fieldInfoDict_.subOrEmptyDict(fieldName);
169
170 scalarMinMax limits;
171
172 fieldDict.readIfPresent("min", limits.min());
173 fieldDict.readIfPresent("max", limits.max());
174
175 return limits;
176}
177
178
180Foam::coordSetWriters::gltfWriter::getAlphaField
181(
182 const dictionary& dict
183) const
184{
185 // Fallback value
186 scalar alphaValue(1);
187
188 const entry* eptr = dict.findEntry("alpha", keyType::LITERAL);
189
190 if (!eptr)
191 {
192 // Not specified
193 }
194 else if (!eptr->stream().peek().isString())
195 {
196 // Value specified
197
198 ITstream& is = eptr->stream();
199 is >> alphaValue;
200 dict.checkITstream(is, "alpha");
201 }
202 else
203 {
204 // Enumeration
205
206 const auto option = fieldOptionNames_.get("alpha", dict);
207
208 switch (option)
209 {
210 case fieldOption::NONE:
211 {
212 break;
213 }
214 case fieldOption::UNIFORM:
215 {
216 dict.readEntry("alphaValue", alphaValue);
217 break;
218 }
219 case fieldOption::FIELD:
220 {
222 << "Unsupported 'field' specification for alpha values"
223 << endl;
224 break;
225 }
226 }
227 }
228
229 return tmp<scalarField>::New(1, alphaValue);
230}
231
232
233void Foam::coordSetWriters::gltfWriter::setupAnimationColour()
234{
235 const dictionary& dict = animationDict_;
236
237 const entry* eptr = dict.findEntry("colour", keyType::LITERAL);
238
239 if (!eptr || !eptr->isStream())
240 {
242 << "Missing 'colour' entry"
243 << exit(FatalIOError);
244 }
245 else if (!eptr->stream().peek().isString())
246 {
247 // Value specified
248
249 ITstream& is = eptr->stream();
250 is >> animateColourValue_;
251 dict.checkITstream(is, "colour");
252
253 // Has uniform value
254 animateColourOption_ = fieldOption::UNIFORM;
255 }
256 else
257 {
258 // Enumeration
259
260 const auto option = fieldOptionNames_.get("colour", dict);
261
262 switch (option)
263 {
264 case fieldOption::NONE:
265 {
267 << "Cannot select 'none' for colour entry!" << nl
268 << "... possible programming error"
269 << exit(FatalError);
270 break;
271 }
272 case fieldOption::UNIFORM:
273 {
274 dict.readEntry("colourValue", animateColourValue_);
275
276 // Has uniform value
277 animateColourOption_ = fieldOption::UNIFORM;
278 break;
279 }
280 case fieldOption::FIELD:
281 {
282 // Needs named field...
283 animateColourName_ = dict.get<word>("colourField");
284 animateColourOption_ = fieldOption::FIELD;
285 break;
286 }
287 }
288 }
289}
290
291
292// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
293
295:
297 writer_(nullptr),
298 animate_(false),
299 colour_(false),
300 animateColourOption_(fieldOption::NONE),
301 animateColourName_(),
302 animateColourValue_(Zero),
303 fieldInfoDict_(),
304 animationDict_()
305{}
306
307
309:
310 coordSetWriter(options),
311 writer_(nullptr),
312 animate_(options.getOrDefault("animate", false)),
313 colour_(options.getOrDefault("colour", false)),
314 animateColourOption_(fieldOption::NONE),
315 animateColourName_(),
316 animateColourValue_(Zero),
317 fieldInfoDict_(options.subOrEmptyDict("fieldInfo")),
318 animationDict_(options.subOrEmptyDict("animationInfo"))
319{
320 // fieldInfo
321 // {
322 // U
323 // {
324 // colourMap coolToWarm; // others...
325 // min 10;
326 // max 100;
327 // alpha 0.5;
328 // }
329 // }
330}
331
332
334(
335 const coordSet& coords,
336 const fileName& outputPath,
337 const dictionary& options
338)
339:
340 gltfWriter(options)
341{
342 open(coords, outputPath);
343}
344
345
347(
348 const UPtrList<coordSet>& tracks,
349 const fileName& outputPath,
350 const dictionary& options
351)
352:
353 gltfWriter(options)
354{
355 open(tracks, outputPath);
356}
357
358
359// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
360
362{
363 close();
364}
365
366
367// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
368
370{
371 // 1) rootdir/<TIME>/setName.gltf
372 // 2) rootdir/setName.gltf
373
374 return getExpectedPath("gltf");
375}
376
377
379{
380 writer_.reset(nullptr);
382}
383
384
386{
387 writer_.reset(nullptr);
389}
390
391
393{
394 writer_.reset(nullptr);
396}
397
398
400{
401 writer_.reset(nullptr);
403}
404
405
406// * * * * * * * * * * * * * * * Implementation * * * * * * * * * * * * * * * //
407
408template<class Type>
409Foam::fileName Foam::coordSetWriters::gltfWriter::writeTemplate
410(
411 const word& fieldName,
412 const UPtrList<const Field<Type>>& fieldPtrs
413)
414{
415 if (coords_.size() != fieldPtrs.size())
416 {
418 << "Attempted to write field: " << fieldName
419 << " (" << fieldPtrs.size() << " entries) for "
420 << coords_.size() << " sets" << nl
421 << exit(FatalError);
422 }
423
424 const auto& tracks = coords_;
425
426 // const auto& times = trackTimes_;
427
428 if (!writer_)
429 {
430 // Field:
431 // 1) rootdir/<TIME>/setName.gltf
432 // 2) rootdir/setName.gltf
433
434 fileName outputFile = path();
435
436 writer_.reset(new glTF::sceneWriter(outputFile));
437
438 auto& scene = writer_->getScene();
439
440 meshes_.resize(tracks.size());
441
442 forAll(tracks, tracki)
443 {
444 word meshName("track:" + Foam::name(tracki));
445 if (tracks.size() == 1)
446 {
447 meshName = "points";
448 }
449
450 meshes_[tracki] = scene.addMesh(tracks[tracki], meshName);
451 }
452 }
453
454
455 auto& scene = writer_->getScene();
456
457 forAll(tracks, tracki)
458 {
459 const label meshi = meshes_[tracki];
460 const auto& field = fieldPtrs[tracki];
461
462 scene.addFieldToMesh(field, fieldName, meshi);
463
464 if (colour_)
465 {
466 const dictionary dict = fieldInfoDict_.subOrEmptyDict(fieldName);
467 const auto& colours = getColourTable(dict);
468
469 const auto talpha = getAlphaField(dict);
470 const scalarField& alpha = talpha();
471
472 const scalarMinMax valLimits = getFieldLimits(fieldName);
473
474 scalarMinMax fldLimits;
476 {
477 MinMax<Type> scanned(minMax(field));
478
479 fldLimits.min() = scalar(component(scanned.min(), 0));
480 fldLimits.max() = scalar(component(scanned.max(), 0));
481 }
482 else
483 {
484 // Use mag() for multiple components
485 fldLimits = minMaxMag(field);
486 }
487
488 // Generated field colours
489 vectorField fieldColour
490 (
492 (
493 colours,
494 field,
495 max(fldLimits.min(), valLimits.min()), // boundMin
496 min(fldLimits.max(), valLimits.max()) // boundMax
497 )
498 );
499
500 scene.addColourToMesh
501 (
502 fieldColour,
503 "Colour:" + fieldName,
504 meshi,
505 alpha
506 );
507 }
508 }
509
510 return writer_().path();
511}
512
513
514template<class Type>
515Foam::fileName Foam::coordSetWriters::gltfWriter::writeTemplate_animate
516(
517 const word& fieldName,
518 const UPtrList<const Field<Type>>& fieldPtrs
519)
520{
521 if (coords_.size() != fieldPtrs.size())
522 {
524 << "Attempted to write field: " << fieldName
525 << " (" << fieldPtrs.size() << " entries) for "
526 << coords_.size() << " sets" << nl
527 << exit(FatalError);
528 }
529
530 const auto& tracks = this->coords_;
531 const auto& times = this->trackTimes_;
532
533 if (!writer_)
534 {
535 // Field:
536 // 1) rootdir/<TIME>/setName.gltf
537 // 2) rootdir/setName.gltf
538
539 fileName outputFile = path();
540
541 writer_.reset(new glTF::sceneWriter(outputFile));
542
543 auto& scene = writer_->getScene();
544
545 meshes_.resize(tracks.size());
546
547 const label animationi = scene.createAnimation("animation");
548
549 forAll(tracks, tracki)
550 {
551 const auto& track = tracks[tracki];
552
553 if (track.empty())
554 {
555 meshes_[tracki] = -1;
556 continue;
557 }
558
559 // Seed starting position
560
561 meshes_[tracki] =
562 scene.addMesh
563 (
564 vectorField(1, track[0]),
565 "track:" + Foam::name(tracki)
566 );
567
568 const label meshi = meshes_[tracki];
569
570 // Time frames
571 const label timeId =
572 scene.addField(times[tracki], "time:" + Foam::name(tracki));
573
574 // Translations
575 const vectorField translation(track - track[0]);
576 const label translationId =
577 scene.addField(translation, "translation");
578
579 scene.addToAnimation(animationi, timeId, translationId, meshi);
580 }
581 }
582
583
584 auto& scene = writer_->getScene();
585
586 // Seed starting field values
587
588 forAll(tracks, tracki)
589 {
590 const auto& track = tracks[tracki];
591 const label meshi = meshes_[tracki];
592 const Field<Type>& field = fieldPtrs[tracki];
593
594 if (track.empty() || meshi < 0)
595 {
596 continue;
597 }
598
599 // Seed starting field values
600 scene.addFieldToMesh(Field<Type>(1, field[0]), fieldName, meshi);
601 }
602
603
604 // Note: colours cannot be animated... setting a fixed value.
605 // However, we need to wait until the field is actually seen
606
607 if (colour_)
608 {
609 if (animateColourOption_ == fieldOption::NONE)
610 {
611 // First time - scan for information
612 setupAnimationColour();
613 }
614
615 switch (animateColourOption_)
616 {
617 case fieldOption::NONE:
618 {
619 // Should not occur
620 break;
621 }
622 case fieldOption::UNIFORM:
623 {
624 // Colour value is known
625
626 vectorField fieldColour(1, animateColourValue_);
627 scalarField alphaChannel(1, 1.0);
628
629 const auto talpha = getAlphaField(animationDict_);
630
631 if (talpha && talpha().size())
632 {
633 alphaChannel[0] = talpha()[0];
634 }
635
636 forAll(tracks, tracki)
637 {
638 const auto& track = tracks[tracki];
639 const label meshi = meshes_[tracki];
640
641 if (track.empty() || meshi < 0)
642 {
643 continue;
644 }
645
646 scene.addColourToMesh
647 (
648 fieldColour,
649 "Colour:fixed", // ... or "Colour:constant"
650 meshi,
651 alphaChannel
652 );
653 }
654
655 // Mark as done
656 animateColourName_.clear();
657 animateColourOption_ = fieldOption::FIELD;
658 break;
659 }
660 case fieldOption::FIELD:
661 {
662 if
663 (
664 !animateColourName_.empty()
665 && animateColourName_ == fieldName
666 )
667 {
668 // This is the desired colour field. Process now
669
670 const auto& colours = getColourTable(animationDict_);
671
672 vectorField fieldColour(1, Zero);
673 scalarField alphaChannel(1, 1.0);
674
675 const auto talpha = getAlphaField(animationDict_);
676
677 if (talpha && talpha().size())
678 {
679 alphaChannel[0] = talpha()[0];
680 }
681
682 forAll(tracks, tracki)
683 {
684 const auto& track = tracks[tracki];
685 const label meshi = meshes_[tracki];
686 const Field<Type>& field = fieldPtrs[tracki];
687
688 if (track.empty() || meshi < 0)
689 {
690 continue;
691 }
692
693 fieldColour[0] =
694 getAnimationColour(animationDict_, colours, field);
695
696 scene.addColourToMesh
697 (
698 fieldColour,
699 "Colour:fixed", // ... or "Colour:constant"
700 meshi,
701 alphaChannel
702 );
703 }
704
705 // Mark colouring as done. Avoid retriggering
706 animateColourName_.clear();
707 animateColourOption_ = fieldOption::FIELD;
708 }
709 break;
710 }
711 }
712 }
713
714 return writer_().path();
715}
716
717
718template<class Type>
719Foam::fileName Foam::coordSetWriters::gltfWriter::writeTemplate
720(
721 const word& fieldName,
722 const Field<Type>& values
723)
724{
725 checkOpen();
726 if (coords_.empty())
727 {
728 return fileName::null;
729 }
730
731 UPtrList<const Field<Type>> fieldPtrs(repackageFields(values));
732 return writeTemplate(fieldName, fieldPtrs);
733}
734
735
736template<class Type>
737Foam::fileName Foam::coordSetWriters::gltfWriter::writeTemplate
738(
739 const word& fieldName,
740 const List<Field<Type>>& fieldValues
741)
742{
743 checkOpen();
744 if (coords_.empty())
745 {
746 return fileName::null;
747 }
748
749 UPtrList<const Field<Type>> fieldPtrs(repackageFields(fieldValues));
750
751 if (animate_ && trackTimes_.size() >= coords_.size())
752 {
753 return writeTemplate_animate(fieldName, fieldPtrs);
754 }
755
756 return writeTemplate(fieldName, fieldPtrs);
757}
758
759
760// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
761
762// Field writing methods
764
765
766// ************************************************************************* //
static float colourMap[]
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
const List< word > & names() const noexcept
The list of enum names, in construction order. Same as toc()
Definition: EnumI.H:46
Generic templated field type.
Definition: Field.H:82
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
const T & max() const noexcept
The max value (second)
Definition: MinMaxI.H:121
const T & min() const noexcept
The min value (first)
Definition: MinMaxI.H:107
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
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
Base class for generating a colour table from node points.
Definition: colourTable.H:80
vector value(const scalar x) const
Return the colour at x. The input is clipped to 0-1 range.
Definition: colourTable.C:94
static const Enum< predefinedType > predefinedNames
Enumeration names for predefinedType.
Definition: colourTable.H:108
Base class for writing coordSet(s) and tracks with fields.
virtual void open(const fileName &outputPath)
Write separate geometry to file.
virtual void endTime()
End a time-step.
virtual void beginTime(const Time &t)
Begin a time-step.
A coordSet(s) writer in glTF v2 format, which is particularly useful for writing track data.
fieldOption
Field option used for colours.
virtual void endTime()
End time step. Clears existing backend.
virtual void beginTime(const Time &t)
Begin time step. Clears existing backend.
static const Enum< fieldOption > fieldOptionNames_
Strings corresponding to the field options.
virtual fileName path() const
Expected (characteristic) output file name - information only.
virtual ~gltfWriter()
Destructor. Calls close()
Holds list of sampling positions.
Definition: coordSet.H:56
reference ref() const
A reference to the entry (Error if not found)
Definition: dictionary.H:225
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
A class for handling file names.
Definition: fileName.H:76
static const fileName null
An empty fileName.
Definition: fileName.H:102
Wrapper for glTF scene for file output.
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
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeName(Type)
Define the typeName.
Definition: className.H:96
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
rDeltaTY field()
Convenience macros for instantiating coordSetWriter methods.
#define defineCoordSetWriterWriteFields(ThisClass)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
static vector getAnimationColour(const dictionary &dict, const colourTable &colours, const Field< Type > &field)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:117
static tmp< vectorField > getBoundedColours(const colourTable &colours, const Field< Type > &field, const scalar boundMin, const scalar boundMax)
MinMax< label > minMax(const labelHashSet &set)
Find the min/max values of labelHashSet.
Definition: hashSets.C:61
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dimensioned< scalarMinMax > minMaxMag(const DimensionedField< Type, GeoMesh > &df)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static constexpr char close
Definition: FlatOutput.H:74