gltfSetWriter.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 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 "gltfSetWriter.H"
29 #include "coordSet.H"
30 #include "fileName.H"
31 #include "OFstream.H"
32 #include "floatVector.H"
33 #include "foamGltfScene.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 template<class Type>
41 ({
42  { fieldOption::UNIFORM, "uniform" },
43  { fieldOption::FIELD, "field" },
44 });
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 template<class Type>
51 (
52  const dictionary& dict
53 ) const
54 {
55  word colourMap = colourTable::predefinedNames.names()[0];
56  dict.readIfPresent("colourMap", colourMap);
57 
58  return colourMap;
59 }
60 
61 
62 template<class Type>
64 (
65  const dictionary& dict
66 ) const
67 {
68  return colourTable::ref(getColourMap(dict));
69 }
70 
71 
72 template<class Type>
74 (
75  const word& fieldName
76 ) const
77 {
78  const dictionary fieldDict = fieldInfoDict_.subOrEmptyDict(fieldName);
79 
80  return fieldDict.getOrDefault("min", -GREAT);
81 }
82 
83 
84 template<class Type>
86 (
87  const word& fieldName
88 ) const
89 {
90  const dictionary fieldDict = fieldInfoDict_.subOrEmptyDict(fieldName);
91 
92  return fieldDict.getOrDefault("max", GREAT);
93 }
94 
95 
96 template<class Type>
98 (
99  const dictionary& dict,
100  const wordList& valueSetNames,
101  const List<const Field<Type>*>& valueSets
102 ) const
103 {
104  if (dict.found("alpha"))
105  {
106  const auto option = fieldOptionNames_.get("alpha", dict);
107 
108  switch (option)
109  {
110  case fieldOption::UNIFORM:
111  {
112  const scalar value = dict.getScalar("alphaValue");
113  return tmp<scalarField>::New(valueSets[0]->size(), value);
114  }
115  case fieldOption::FIELD:
116  {
117  const word alphaFieldName = dict.get<word>("alphaField");
118  const bool normalise = dict.get<bool>("normalise");
119  const label i = valueSetNames.find(alphaFieldName);
120  if (i == -1)
121  {
123  << "Unable to find field " << alphaFieldName
124  << ". Valid field names are:" << valueSetNames
125  << exit(FatalError);
126  }
127 
128  auto tresult =
129  tmp<scalarField>::New(valueSets[i]->component(0));
130 
131  if (normalise)
132  {
133  tresult.ref() /= mag(tresult() + ROOTVSMALL);
134  }
135 
136  return tresult;
137  }
138  }
139  }
140 
141  return tmp<scalarField>::New(valueSets[0]->size(), Zero);
142 }
143 
144 
145 template<class Type>
147 (
148  const dictionary& dict,
149  const wordList& valueSetNames,
150  const List<List<Field<Type>>>& valueSets,
151  const label tracki
152 ) const
153 {
154  if (dict.found("alpha"))
155  {
156  const auto option = fieldOptionNames_.get("alpha", dict);
157 
158  switch (option)
159  {
160  case fieldOption::UNIFORM:
161  {
162  const scalar value = dict.getScalar("alphaValue");
163  return tmp<scalarField>::New
164  (
165  valueSets[0][tracki].size(), value
166  );
167  }
168  case fieldOption::FIELD:
169  {
170  const word alphaFieldName = dict.get<word>("alphaField");
171  const bool normalise = dict.get<bool>("normalise");
172  const label fieldi = valueSetNames.find(alphaFieldName);
173  if (fieldi == -1)
174  {
176  << "Unable to find field " << alphaFieldName
177  << ". Valid field names are:" << valueSetNames
178  << exit(FatalError);
179  }
180 
181  // Note: selecting the first component!
182  auto tresult =
184  (
185  valueSets[fieldi][tracki].component(0)
186  );
187 
188  if (normalise)
189  {
190  tresult.ref() /= mag(tresult() + ROOTVSMALL);
191  }
192 
193  return tresult;
194  }
195  }
196  }
197 
198  return tmp<scalarField>::New(valueSets[0][tracki].size(), Zero);
199 }
200 
201 
202 template<class Type>
204 (
205  const colourTable& colours,
206  const wordList& valueSetNames,
207  const List<List<Field<Type>>>& valueSets,
208  const label tracki
209 ) const
210 {
211  if (!colour_)
212  {
214  << "Attempting to get colour when colour option is off"
215  << abort(FatalError);
216  }
217 
218  const auto option = fieldOptionNames_.get("colour", animationDict_);
219 
220  switch (option)
221  {
222  case fieldOption::UNIFORM:
223  {
224  return animationDict_.get<vector>("colourValue");
225  }
226  case fieldOption::FIELD:
227  {
228  const word fieldName = animationDict_.get<word>("colourField");
229  const label fieldi = valueSetNames.find(fieldName);
230  if (fieldi == -1)
231  {
233  << "Unable to find field " << fieldName
234  << ". Valid field names are:" << valueSetNames
235  << exit(FatalError);
236  }
237 
238  // Note: selecting the first component!
239 
240  scalar minValue;
241  scalar maxValue;
242  if (!animationDict_.readIfPresent("min", minValue))
243  {
244  minValue = min(valueSets[fieldi][tracki].component(0));
245  }
246  if (!animationDict_.readIfPresent("max", maxValue))
247  {
248  maxValue = max(valueSets[fieldi][tracki].component(0));
249  }
250  const scalar refValue = component(valueSets[fieldi][tracki][0], 0);
251  const scalar fraction =
252  (refValue - minValue)/(maxValue - minValue + ROOTVSMALL);
253 
254  return (colours.value(max(0, min(1, fraction))));
255  }
256  }
257 
258  return vector::zero;
259 }
260 
261 
262 template<class Type>
264 (
265  const coordSet& points
266 ) const
267 {
268  auto tresult = tmp<vectorField>::New(points.size(), Zero);
269  auto& result = tresult.ref();
270 
271  if (points.size() > 1)
272  {
273  for (label i = 1; i < points.size(); ++i)
274  {
275  result[i-1] = points[i] - points[i-1];
276  result[i-1].normalise();
277  }
278 
279  result.last() = result[points.size()-2];
280  }
281 
282 
283  return tresult;
284 }
285 
286 
287 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
288 
289 template<class Type>
291 :
292  writer<Type>(),
293  animate_(false),
294  colour_(false),
295  fieldInfoDict_(),
296  animationDict_()
297 {}
298 
299 
300 template<class Type>
302 :
303  writer<Type>(dict),
304  animate_(dict.getOrDefault("animate", false)),
305  colour_(dict.getOrDefault("colour", false)),
306  fieldInfoDict_(dict.subOrEmptyDict("fieldInfo")),
307  animationDict_(dict.subOrEmptyDict("animationInfo"))
308 {
309  // fieldInfo
310  // {
311  // U
312  // {
313  // colourMap coolToWarm; // others...
314  // min 10;
315  // max 100;
316  // alpha field; // uniform|field
317  // alphaField ageOfAir;
318  // }
319  // }
320 }
321 
322 
323 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
324 
325 template<class Type>
327 (
328  const coordSet& points,
329  const wordList& valueSetNames
330 ) const
331 {
332  return this->getBaseName(points, valueSetNames) + ".gltf";
333 }
334 
335 
336 template<class Type>
338 (
339  const coordSet& points,
340  const wordList& valueSetNames,
341  const List<const Field<Type>*>& valueSets,
342  Ostream& os
343 ) const
344 {
345  if (valueSets.size() != valueSetNames.size())
346  {
348  << "Number of variables:" << valueSetNames.size() << endl
349  << "Number of valueSets:" << valueSets.size()
350  << exit(FatalError);
351  }
352 
353  glTF::scene scene;
354  const label meshi = scene.addMesh(points, "points");
355  forAll(valueSetNames, i)
356  {
357  scene.addFieldToMesh(*valueSets[i], valueSetNames[i], meshi);
358  }
359 
360  if (colour_)
361  {
362  forAll(valueSets, fieldi)
363  {
364  const auto& field = *valueSets[fieldi];
365  const word& fieldName = valueSetNames[fieldi];
366  const dictionary dict = fieldInfoDict_.subOrEmptyDict(fieldName);
367  const auto& colours = getColourTable(dict);
368 
369  const auto talpha =
370  getAlphaField(dict, valueSetNames, valueSets);
371  const scalarField& alpha = talpha();
372 
373  const Type maxValue = max(field);
374  const Type minValue = min(field);
375 
376  const scalar minValueLimit = getFieldMin(fieldName);
377  const scalar maxValueLimit = getFieldMax(fieldName);
378 
379  for (direction cmpti=0; cmpti < pTraits<Type>::nComponents; ++cmpti)
380  {
381  vectorField fieldColour(field.size());
382 
383  forAll(field, i)
384  {
385  const Type& v = field[i];
386  float f = component(v, cmpti);
387  float minf = max(component(minValue, cmpti), minValueLimit);
388  float maxf = min(component(maxValue, cmpti), maxValueLimit);
389  float deltaf = (maxf - minf + SMALL);
390 
391  fieldColour[i] =
392  colours.value(min(max((f - minf)/deltaf, 0), 1));
393  }
394 
395  scene.addColourToMesh
396  (
397  fieldColour,
398  "Colour:" + fieldName + Foam::name(cmpti),
399  meshi,
400  alpha
401  );
402  }
403  }
404  }
405 
406  scene.write(os);
407 }
408 
409 
410 template<class Type>
412 (
413  const bool writeTracks,
414  const List<scalarField>& times,
415  const PtrList<coordSet>& tracks,
416  const wordList& valueSetNames,
417  const List<List<Field<Type>>>& valueSets,
418  Ostream& os
419 ) const
420 {
421  if (valueSets.size() != valueSetNames.size())
422  {
424  << "Number of variables:" << valueSetNames.size() << endl
425  << "Number of valueSets:" << valueSets.size()
426  << exit(FatalError);
427  }
428 
429  if (animate_)
430  {
431  writeAnimateTracks
432  (
433  writeTracks,
434  times,
435  tracks,
436  valueSetNames,
437  valueSets,
438  os
439  );
440  }
441  else
442  {
443  writeStaticTracks
444  (
445  writeTracks,
446  times,
447  tracks,
448  valueSetNames,
449  valueSets,
450  os
451  );
452  }
453 }
454 
455 
456 template<class Type>
458 (
459  const bool writeTracks,
460  const List<scalarField>& times,
461  const PtrList<coordSet>& tracks,
462  const wordList& valueSetNames,
463  const List<List<Field<Type>>>& valueSets,
464  Ostream& os
465 ) const
466 {
467  glTF::scene scene;
468  forAll(tracks, tracki)
469  {
470  const vectorField& track = tracks[tracki];
471  const label meshi = scene.addMesh(track, "track:" + Foam::name(tracki));
472  forAll(valueSetNames, fieldi)
473  {
474  const word& fieldName = valueSetNames[fieldi];
475  const auto& field = valueSets[fieldi][tracki];
476  scene.addFieldToMesh(field, fieldName, meshi);
477  }
478 
479  if (colour_)
480  {
481  forAll(valueSets, fieldi)
482  {
483  const auto& field = valueSets[fieldi][tracki];
484  const word& fieldName = valueSetNames[fieldi];
485  const dictionary dict =
486  fieldInfoDict_.subOrEmptyDict(fieldName);
487  const auto& colours = getColourTable(dict);
488 
489  const auto talpha =
490  getTrackAlphaField(dict, valueSetNames, valueSets, tracki);
491  const scalarField& alpha = talpha();
492 
493  const Type maxValue = max(field);
494  const Type minValue = min(field);
495 
496  const scalar minValueLimit = getFieldMin(fieldName);
497  const scalar maxValueLimit = getFieldMax(fieldName);
498 
499  for
500  (
501  direction cmpti=0;
502  cmpti < pTraits<Type>::nComponents;
503  ++cmpti
504  )
505  {
506  vectorField fieldColour(field.size(), Zero);
507 
508  forAll(field, i)
509  {
510  const Type& v = field[i];
511  float f = component(v, cmpti);
512  float minf =
513  max(component(minValue, cmpti), minValueLimit);
514  float maxf =
515  min(component(maxValue, cmpti), maxValueLimit);
516  float deltaf = (maxf - minf + SMALL);
517 
518  fieldColour[i] =
519  colours.value(min(max((f - minf)/deltaf, 0), 1));
520  }
521 
522  scene.addColourToMesh
523  (
524  fieldColour,
525  "Colour:" + fieldName + Foam::name(cmpti),
526  meshi,
527  alpha
528  );
529  }
530  }
531  }
532  }
533 
534  scene.write(os);
535 }
536 
537 
538 template<class Type>
540 (
541  const bool writeTracks,
542  const List<scalarField>& times,
543  const PtrList<coordSet>& tracks,
544  const wordList& valueSetNames,
545  const List<List<Field<Type>>>& valueSets,
546  Ostream& os
547 ) const
548 {
549  const auto& colours = getColourTable(animationDict_);
550 
551  glTF::scene scene;
552  const label animationi = scene.createAnimation("animation");
553 
554  forAll(tracks, tracki)
555  {
556  const auto& track = tracks[tracki];
557 
558  if (track.empty())
559  {
560  continue;
561  }
562 
563  // Seed starting positions and field values
564  const label meshi =
565  scene.addMesh
566  (
567  vectorField(1, track[0]),
568  "track:" + Foam::name(tracki)
569  );
570 
571  forAll(valueSetNames, fieldi)
572  {
573  const Field<Type>& field = valueSets[fieldi][tracki];
574  const word& fieldName = valueSetNames[fieldi];
575  scene.addFieldToMesh(Field<Type>(1, field[0]), fieldName, meshi);
576  }
577 
578  // Time frames
579  const label timeId =
580  scene.addField(times[tracki], "time:" + Foam::name(tracki));
581 
582  // Translations
583  const vectorField translation(track - track[0]);
584  const label translationId = scene.addField(translation, "translation");
585 
586  scene.addToAnimation(animationi, timeId, translationId, meshi);
587 
588  // Note: colours cannot be animated... setting a fixed value
589  if (colour_)
590  {
591  const vector colour =
592  getTrackAnimationColour
593  (
594  colours,
595  valueSetNames,
596  valueSets,
597  tracki
598  );
599 
600  const auto talpha =
601  getTrackAlphaField
602  (
603  animationDict_,
604  valueSetNames,
605  valueSets,
606  tracki
607  );
608 
609  const scalarField& alpha = talpha();
610 
611  scene.addColourToMesh
612  (
613  vectorField(1, colour),
614  "Colour:fixed",
615  meshi,
616  scalarField(1, alpha[0])
617  );
618  }
619  }
620 
621  scene.write(os);
622 }
623 
624 
625 // ************************************************************************* //
floatVector.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:44
Foam::gltfSetWriter::gltfSetWriter
gltfSetWriter()
Default construct.
Definition: gltfSetWriter.C:290
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::gltfSetWriter::fieldOptionNames_
static const Enum< fieldOption > fieldOptionNames_
Strings corresponding to the field options.
Definition: gltfSetWriter.H:161
coordSet.H
Foam::glTF::scene::addField
label addField(const Type &fld, const word &name, const label target=-1)
Returns accessor index.
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::glTF::scene::addMesh
label addMesh(const Type &fld, const word &name)
Returns index of last mesh.
Foam::gltfSetWriter::getFileName
virtual fileName getFileName(const coordSet &, const wordList &) const
Return the file name.
Definition: gltfSetWriter.C:327
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
foamGltfScene.H
Foam::glTF::scene::addToAnimation
void addToAnimation(const label animationi, const label inputId, const label outputId, const label meshId, const string &interpolation="LINEAR")
Add to existing animation.
Definition: foamGltfScene.C:102
ref
rDeltaT ref()
Foam::gltfSetWriter
Writes point data in glTF v2 format.
Definition: gltfSetWriter.H:144
gltfSetWriter.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::gltfSetWriter::writeStaticTracks
virtual void writeStaticTracks(const bool writeTracks, const List< scalarField > &times, const PtrList< coordSet > &tracks, const wordList &valueSetNames, const List< List< Field< Type >>> &valueSets, Ostream &) const
Write static tracks.
Definition: gltfSetWriter.C:458
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
Foam::glTF::scene::addFieldToMesh
label addFieldToMesh(const Type &fld, const word &name, const label meshi)
Returns accessor index.
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::glTF::scene
Main class to assemble glTF components into a scene.
Definition: foamGltfScene.H:64
colourMap
static float colourMap[]
Definition: AC3DsurfaceFormatCore.C:36
fileName.H
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:81
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
Foam::glTF::scene::write
void write(Ostream &os)
Write to stream (JSON and binary data)
Definition: foamGltfScene.C:139
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::gltfSetWriter::writeAnimateTracks
virtual void writeAnimateTracks(const bool writeTracks, const List< scalarField > &times, const PtrList< coordSet > &tracks, const wordList &valueSetNames, const List< List< Field< Type >>> &valueSets, Ostream &) const
Write animated tracks.
Definition: gltfSetWriter.C:540
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< word >
Foam::glTF::scene::addColourToMesh
label addColourToMesh(const vectorField &fld, const word &name, const label meshi, const scalarField &alpha=scalarField())
Returns accessor index.
Definition: foamGltfScene.C:47
Foam::glTF::scene::createAnimation
label createAnimation(const word &name)
Returns index of last animation.
Definition: foamGltfScene.C:94
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::colourTable
Base class for generating a colour table from node points.
Definition: colourTable.H:79
minValue
scalar minValue
Definition: LISASMDCalcMethod2.H:12
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::gltfSetWriter::write
virtual void write(const coordSet &, const wordList &, const List< const Field< Type > * > &, Ostream &) const
Write.
Definition: gltfSetWriter.C:338
maxValue
scalar maxValue
Definition: LISASMDCalcMethod1.H:5