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