vtkPVFoamUpdateInfo.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "vtkPVFoam.H"
27 
28 // OpenFOAM includes
29 #include "vtkPVFoamReader.h"
30 #include "cellSet.H"
31 #include "faceSet.H"
32 #include "pointSet.H"
33 #include "IOobjectList.H"
34 #include "IOPtrList.H"
36 #include "entry.H"
37 #include "Cloud.H"
38 #include "surfaceFields.H"
39 
40 // Local includes
42 
43 // VTK includes
44 #include "vtkDataArraySelection.h"
45 
46 
47 // * * * * * * * * * * * * * * * Private Classes * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 //- A class for reading zone information without requiring a mesh
53 template<class ZonesType>
54 class zonesEntries
55 :
56  public regIOobject,
57  public PtrList<entry>
58 {
59 
60 public:
61 
62  // Constructors
63 
64  explicit zonesEntries(const IOobject& io)
65  :
66  regIOobject(io),
67  PtrList<entry>(readStream(ZonesType::typeName))
68  {
69  close();
70  }
71 
72  // Member functions
73 
74  bool writeData(Ostream&) const
75  {
77  return false;
78  }
79 };
80 
81 }
82 
83 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
84 
85 template<class ZonesType>
86 Foam::wordList Foam::vtkPVFoam::getZoneNames(const ZonesType& zones) const
87 {
88  wordList names(zones.size());
89  label nZone = 0;
90 
91  forAll(zones, zoneI)
92  {
93  if (zones[zoneI].size())
94  {
95  names[nZone++] = zones[zoneI].name();
96  }
97  }
98  names.setSize(nZone);
99 
100  return names;
101 }
102 
103 
104 template<class ZonesType>
105 Foam::wordList Foam::vtkPVFoam::getZoneNames(const word& zonesName) const
106 {
107  wordList names;
108 
109  // Mesh not loaded - read from file
110  typeIOobject<ZonesType> ioObj
111  (
112  zonesName,
113  dbPtr_().findInstance
114  (
115  meshDir_,
116  zonesName,
118  ),
119  meshDir_,
120  dbPtr_(),
123  false
124  );
125 
126  if (ioObj.headerOk())
127  {
128  const zonesEntries<ZonesType> zones(ioObj);
129 
130  names.setSize(zones.size());
131  forAll(zones, zoneI)
132  {
133  names[zoneI] = zones[zoneI].keyword();
134  }
135  }
136 
137  return names;
138 }
139 
140 
141 void Foam::vtkPVFoam::updateInfoInternalMesh
142 (
143  vtkDataArraySelection* arraySelection
144 )
145 {
146  if (debug)
147  {
148  InfoInFunction << endl;
149  }
150 
151  // Determine mesh parts (internalMesh, patches...)
152  //- Add internal mesh as first entry
153  arrayRangeVolume_.reset(arraySelection->GetNumberOfArrays());
154  arraySelection->AddArray
155  (
156  "internalMesh"
157  );
158  arrayRangeVolume_ += 1;
159 
160  if (debug)
161  {
162  // Just for debug info
163  getSelectedArrayEntries(arraySelection);
164  }
165 }
166 
167 
168 void Foam::vtkPVFoam::updateInfoLagrangian
169 (
170  vtkDataArraySelection* arraySelection
171 )
172 {
173  if (debug)
174  {
175  InfoInFunction << nl
176  << " " << dbPtr_->timePath()/cloud::prefix << endl;
177  }
178 
179 
180  // Use the db directly since this might be called without a mesh,
181  // but the region must get added back in
182  fileName lagrangianPrefix(cloud::prefix);
183  if (meshRegion_ != polyMesh::defaultRegion)
184  {
185  lagrangianPrefix = meshRegion_/cloud::prefix;
186  }
187 
188  arrayRangeLagrangian_.reset(arraySelection->GetNumberOfArrays());
189 
190  // Generate a list of lagrangian clouds across all times
191  HashSet<fileName> cloudDirs;
192 
193  // Get times list. Flush first to force refresh.
194  fileHandler().flush();
195  instantList times = dbPtr_().times();
196  forAll(times, timei)
197  {
198  cloudDirs +=
200  (
201  dbPtr_->path()/times[timei].name()/lagrangianPrefix,
203  );
204  }
205 
206  forAllConstIter(HashSet<fileName>, cloudDirs, cloudIter)
207  {
208  arraySelection->AddArray
209  (
210  (cloudIter.key() + " - lagrangian").c_str()
211  );
212  }
213 
214  arrayRangeLagrangian_ += cloudDirs.size();
215 
216  if (debug)
217  {
218  // Just for debug info
219  getSelectedArrayEntries(arraySelection);
220  }
221 }
222 
223 
224 void Foam::vtkPVFoam::updateInfoPatches
225 (
226  vtkDataArraySelection* arraySelection,
227  stringList& enabledEntries,
228  const bool first
229 )
230 {
231  if (debug)
232  {
234  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "]" << endl;
235  }
236 
237 
238  HashSet<string> enabledEntriesSet(enabledEntries);
239 
240  arrayRangePatches_.reset(arraySelection->GetNumberOfArrays());
241 
242  int nPatches = 0;
243  if (meshPtr_)
244  {
245  const polyBoundaryMesh& patches = meshPtr_->boundaryMesh();
246  const HashTable<labelList, word>& groups = patches.groupPatchIndices();
247  const wordList allPatchNames = patches.names();
248 
249  // Add patch groups
250 
251  for
252  (
253  HashTable<labelList, word>::const_iterator iter = groups.begin();
254  iter != groups.end();
255  ++iter
256  )
257  {
258  const word& groupName = iter.key();
259  const labelList& patchIDs = iter();
260 
261  label nFaces = 0;
262  forAll(patchIDs, i)
263  {
264  nFaces += patches[patchIDs[i]].size();
265  }
266 
267  // Valid patch if nFace > 0 - add patch to GUI list
268  if (nFaces)
269  {
270  string vtkGrpName = groupName + " - group";
271  arraySelection->AddArray(vtkGrpName.c_str());
272 
273  ++nPatches;
274 
275  if (enabledEntriesSet.found(vtkGrpName))
276  {
277  if (!reader_->GetShowGroupsOnly())
278  {
279  enabledEntriesSet.erase(vtkGrpName);
280  forAll(patchIDs, i)
281  {
282  const polyPatch& pp = patches[patchIDs[i]];
283  if (pp.size())
284  {
285  string vtkPatchName
286  (
287  pp.name() + " - " + pp.type()
288  );
289 
290  enabledEntriesSet.insert(vtkPatchName);
291  }
292  }
293  }
294  }
295  }
296  }
297 
298 
299  // Add patches
300 
301  if (!reader_->GetShowGroupsOnly())
302  {
304  {
305  const polyPatch& pp = patches[patchi];
306 
307  if (pp.size())
308  {
309  const string vtkPatchName = pp.name() + " - " + pp.type();
310 
311  // Add patch to GUI list
312  arraySelection->AddArray(vtkPatchName.c_str());
313 
314  ++nPatches;
315  }
316  }
317  }
318  }
319  else
320  {
321  // Mesh not loaded - read from file
322  // but this could fail if we've supplied a bad region name
323  typeIOobject<polyBoundaryMesh> ioObj
324  (
325  "boundary",
326  dbPtr_().findInstance
327  (
328  meshDir_,
329  "boundary",
331  ),
332  meshDir_,
333  dbPtr_(),
336  false
337  );
338 
339  // This should only ever fail if the mesh region doesn't exist
340  if (ioObj.headerOk())
341  {
342  polyBoundaryMeshEntries patchEntries(ioObj);
343 
344 
345  // Read patches and determine sizes
346 
347  wordList names(patchEntries.size());
348  labelList sizes(patchEntries.size());
349 
350  forAll(patchEntries, patchi)
351  {
352  const dictionary& patchDict = patchEntries[patchi].dict();
353 
354  sizes[patchi] = patchDict.lookup<label>("nFaces");
355  names[patchi] = patchEntries[patchi].keyword();
356  }
357 
358 
359  // Add (non-zero) patch groups to the list of mesh parts
360 
361  HashTable<labelList, word> groups(patchEntries.size());
362 
363  forAll(patchEntries, patchi)
364  {
365  const dictionary& patchDict = patchEntries[patchi].dict();
366 
367  wordList groupNames;
368  patchDict.readIfPresent("inGroups", groupNames);
369 
370  forAll(groupNames, groupI)
371  {
372  HashTable<labelList, word>::iterator iter = groups.find
373  (
374  groupNames[groupI]
375  );
376  if (iter != groups.end())
377  {
378  iter().append(patchi);
379  }
380  else
381  {
382  groups.insert(groupNames[groupI], labelList(1, patchi));
383  }
384  }
385  }
386 
387  for
388  (
390  groups.begin();
391  iter != groups.end();
392  ++iter
393  )
394  {
395  const word& groupName = iter.key();
396  const labelList& patchIDs = iter();
397 
398  label nFaces = 0;
399  forAll(patchIDs, i)
400  {
401  nFaces += sizes[patchIDs[i]];
402  }
403 
404  // Valid patch if nFace > 0 - add patch to GUI list
405  if (nFaces)
406  {
407  string vtkGrpName = groupName + " - group";
408  arraySelection->AddArray(vtkGrpName.c_str());
409 
410  ++nPatches;
411 
412  if (enabledEntriesSet.found(vtkGrpName))
413  {
414  if (!reader_->GetShowGroupsOnly())
415  {
416  enabledEntriesSet.erase(vtkGrpName);
417  forAll(patchIDs, i)
418  {
419  if (sizes[patchIDs[i]])
420  {
421  const word patchType
422  (
423  patchEntries[patchIDs[i]].dict().lookup
424  (
425  "type"
426  )
427  );
428 
429  string vtkPatchName
430  (
431  names[patchIDs[i]] + " - " + patchType
432  );
433 
434  enabledEntriesSet.insert(vtkPatchName);
435  }
436  }
437  }
438  }
439  }
440  }
441 
442 
443  // Add (non-zero) patches to the list of mesh parts
444 
445  if (!reader_->GetShowGroupsOnly())
446  {
447  const wordReList defaultPatchTypes
448  (
449  configDict_.lookupOrDefault
450  (
451  "defaultPatchTypes",
452  wordReList(wordList({"patch", "wall"}))
453  )
454  );
455 
456  forAll(names, patchi)
457  {
458  // Valid patch if nFace > 0 - add patch to GUI list
459  if (sizes[patchi])
460  {
461  const word patchType
462  (
463  patchEntries[patchi].dict().lookup("type")
464  );
465 
466  const string vtkPatchName
467  (
468  names[patchi] + " - " + patchType
469  );
470 
471  arraySelection->AddArray(vtkPatchName.c_str());
472 
473  if (first)
474  {
475  if (findStrings(defaultPatchTypes, patchType))
476  {
477  enabledEntriesSet.insert(vtkPatchName);
478  }
479  }
480 
481  ++nPatches;
482  }
483  }
484  }
485  }
486  }
487  arrayRangePatches_ += nPatches;
488 
489  // Update enabled entries in case of group selection
490  enabledEntries = enabledEntriesSet.toc();
491 
492  if (debug)
493  {
494  // Just for debug info
495  getSelectedArrayEntries(arraySelection);
496  }
497 }
498 
499 
500 void Foam::vtkPVFoam::updateInfoZones
501 (
502  vtkDataArraySelection* arraySelection
503 )
504 {
505  if (!reader_->GetIncludeZones())
506  {
507  return;
508  }
509 
510  if (debug)
511  {
513  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "]" << endl;
514  }
515 
516  wordList namesLst;
517 
518  // cellZones information
519  if (meshPtr_)
520  {
521  namesLst = getZoneNames(meshPtr_->cellZones());
522  }
523  else
524  {
525  namesLst = getZoneNames<cellZoneList>("cellZones");
526  }
527 
528  arrayRangeCellZones_.reset(arraySelection->GetNumberOfArrays());
529  forAll(namesLst, elemI)
530  {
531  arraySelection->AddArray
532  (
533  (namesLst[elemI] + " - cellZone").c_str()
534  );
535  }
536  arrayRangeCellZones_ += namesLst.size();
537 
538 
539  // faceZones information
540  if (meshPtr_)
541  {
542  namesLst = getZoneNames(meshPtr_->faceZones());
543  }
544  else
545  {
546  namesLst = getZoneNames<faceZoneList>("faceZones");
547  }
548 
549  arrayRangeFaceZones_.reset(arraySelection->GetNumberOfArrays());
550  forAll(namesLst, elemI)
551  {
552  arraySelection->AddArray
553  (
554  (namesLst[elemI] + " - faceZone").c_str()
555  );
556  }
557  arrayRangeFaceZones_ += namesLst.size();
558 
559 
560  // pointZones information
561  if (meshPtr_)
562  {
563  namesLst = getZoneNames(meshPtr_->pointZones());
564  }
565  else
566  {
567  namesLst = getZoneNames<pointZoneList>("pointZones");
568  }
569 
570  arrayRangePointZones_.reset(arraySelection->GetNumberOfArrays());
571  forAll(namesLst, elemI)
572  {
573  arraySelection->AddArray
574  (
575  (namesLst[elemI] + " - pointZone").c_str()
576  );
577  }
578  arrayRangePointZones_ += namesLst.size();
579 
580  if (debug)
581  {
582  // Just for debug info
583  getSelectedArrayEntries(arraySelection);
584  }
585 }
586 
587 
588 void Foam::vtkPVFoam::updateInfoSets
589 (
590  vtkDataArraySelection* arraySelection
591 )
592 {
593  if (!reader_->GetIncludeSets())
594  {
595  return;
596  }
597 
598  if (debug)
599  {
600  InfoInFunction << endl;
601  }
602 
603  // Add names of sets. Search for last time directory with a sets
604  // subdirectory. Take care not to search beyond the last mesh.
605 
606  word facesInstance = dbPtr_().findInstance
607  (
608  meshDir_,
609  "faces",
611  );
612 
613  word setsInstance = dbPtr_().findInstance
614  (
615  meshDir_/"sets",
616  word::null,
618  facesInstance
619  );
620 
621  IOobjectList objects(dbPtr_(), setsInstance, meshDir_/"sets");
622 
623  if (debug)
624  {
625  Info<< " Foam::vtkPVFoam::updateInfoSets read "
626  << objects.names() << " from " << setsInstance << endl;
627  }
628 
629 
630  arrayRangeCellSets_.reset(arraySelection->GetNumberOfArrays());
631  arrayRangeCellSets_ += addToSelection<cellSet>
632  (
633  arraySelection,
634  objects,
635  " - cellSet"
636  );
637 
638  arrayRangeFaceSets_.reset(arraySelection->GetNumberOfArrays());
639  arrayRangeFaceSets_ += addToSelection<faceSet>
640  (
641  arraySelection,
642  objects,
643  " - faceSet"
644  );
645 
646  arrayRangePointSets_.reset(arraySelection->GetNumberOfArrays());
647  arrayRangePointSets_ += addToSelection<pointSet>
648  (
649  arraySelection,
650  objects,
651  " - pointSet"
652  );
653 
654  if (debug)
655  {
656  // Just for debug info
657  getSelectedArrayEntries(arraySelection);
658  }
659 }
660 
661 
662 void Foam::vtkPVFoam::updateInfoFields()
663 {
664  if (debug)
665  {
667  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "]"
668  << endl;
669  }
670 
671  vtkDataArraySelection* fieldSelection = reader_->GetFieldSelection();
672 
673  // Use the db directly since this might be called without a mesh,
674  // but the region must get added back in
675  word regionPrefix;
676  if (meshRegion_ != polyMesh::defaultRegion)
677  {
678  regionPrefix = meshRegion_;
679  }
680 
681  const Time& runTime = dbPtr_();
682  const instantList times = runTime.times();
683 
684  stringList enabledEntries;
685 
686  if (fieldSelection->GetNumberOfArrays() == 0 && !meshPtr_)
687  {
688  const wordReList defaultFieldRes
689  (
690  configDict_.lookupOrDefault
691  (
692  "defaultFields",
693  wordReList(wordList{"p", "U"})
694  )
695  );
696 
697  wordHashSet objectNameSet;
698  forAll(times, timei)
699  {
700  // Search for list of objects for this time and mesh region
701  IOobjectList objects(runTime, times[timei].name(), regionPrefix);
702 
703  forAllConstIter(IOobjectList, objects, iter)
704  {
705  objectNameSet.insert(iter.key());
706  }
707  }
708 
709  const wordList objectNames(objectNameSet.toc());
710  const labelList defaultFields
711  (
712  findStrings(defaultFieldRes, objectNames)
713  );
714 
715  enabledEntries.setSize(defaultFields.size());
716  forAll(defaultFields, i)
717  {
718  enabledEntries[i] = objectNames[defaultFields[i]];
719  }
720  }
721  else
722  {
723  // Preserve the enabled selections
724  enabledEntries = getSelectedArrayEntries(fieldSelection);
725  }
726 
727  fieldSelection->RemoveAllArrays();
728 
729  forAll(times, timei)
730  {
731  // Search for list of objects for this time and mesh region
732  IOobjectList objects(runTime, times[timei].name(), regionPrefix);
733 
734  addFieldsToSelection<fvPatchField, volMesh>
735  (
736  fieldSelection,
737  objects
738  );
739 
740  addInternalFieldsToSelection<fvPatchField, volMesh>
741  (
742  fieldSelection,
743  objects
744  );
745 
746  addFieldsToSelection<fvsPatchField, surfaceMesh>
747  (
748  fieldSelection,
749  objects
750  );
751 
752  addFieldsToSelection<pointPatchField, pointMesh>
753  (
754  fieldSelection,
755  objects
756  );
757  }
758 
759  // Restore the enabled selections
760  setSelectedArrayEntries(fieldSelection, enabledEntries);
761 }
762 
763 
764 void Foam::vtkPVFoam::updateInfoLagrangianFields()
765 {
766  if (debug)
767  {
768  InfoInFunction << endl;
769  }
770 
771  vtkDataArraySelection* fieldSelection =
772  reader_->GetLagrangianFieldSelection();
773 
774  // Preserve the enabled selections
775  stringList enabledEntries = getSelectedArrayEntries(fieldSelection);
776  fieldSelection->RemoveAllArrays();
777 
778  // Use the db directly since this might be called without a mesh,
779  // but the region must get added back in
780  fileName lagrangianPrefix(cloud::prefix);
781  if (meshRegion_ != polyMesh::defaultRegion)
782  {
783  lagrangianPrefix = meshRegion_/cloud::prefix;
784  }
785 
786  // Add the available fields from all clouds and all time directories.
787  // Differing sets of fields from multiple clouds get combined into a single
788  // set. ParaView will display "(partial)" after field names that only apply
789  // to some of the clouds.
790  const arrayRange& range = arrayRangeLagrangian_;
791 
792  fileHandler().flush();
793  for (label partId = range.start(); partId < range.end(); ++ partId)
794  {
795  const instantList times = dbPtr_().times();
796  forAll(times, timei)
797  {
798  IOobjectList objects
799  (
800  dbPtr_(),
801  times[timei].name(),
802  lagrangianPrefix/getPartName(partId)
803  );
804 
805  addToSelection<IOField<label>>(fieldSelection, objects);
806  addToSelection<IOField<scalar>>(fieldSelection, objects);
807  addToSelection<IOField<vector>>(fieldSelection, objects);
808  addToSelection<IOField<sphericalTensor>>(fieldSelection, objects);
809  addToSelection<IOField<symmTensor>>(fieldSelection, objects);
810  addToSelection<IOField<tensor>>(fieldSelection, objects);
811  }
812  }
813 
814  // Restore the enabled selections
815  setSelectedArrayEntries(fieldSelection, enabledEntries);
816 }
817 
818 
819 // ************************************************************************* //
scalar range
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
friend class iterator
Declare friendship with the iterator.
Definition: HashTable.H:194
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:197
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:63
virtual void flush() const
Forcibly wait until all output done. Flush any cached data.
virtual fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true) const =0
Read a directory and return the entries as a string list.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:267
static const word null
An empty word.
Definition: word.H:77
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
label patchi
const fvPatchList & patches
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
List< instant > instantList
List of instants.
Definition: instantList.H:42
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
List< string > stringList
A List of strings.
Definition: stringList.H:50
static const char nl
Definition: Ostream.H:266
label nPatches
Definition: readKivaGrid.H:396
objects
dictionary dict
Foam::surfaceFields.
const bool writeData(readBool(pdfDictionary.lookup("writeData")))