vtkPVFoamUpdateInfo.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 | Copyright (C) 2011-2017 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 "cellSet.H"
30 #include "faceSet.H"
31 #include "pointSet.H"
32 #include "IOobjectList.H"
33 #include "IOPtrList.H"
35 #include "entry.H"
36 #include "Cloud.H"
37 #include "vtkPVFoamReader.h"
38 
39 // local headers
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  Info<< "<beg> Foam::vtkPVFoam::updateInfoInternalMesh" << 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  Info<< "<end> Foam::vtkPVFoam::updateInfoInternalMesh" << endl;
167  }
168 }
169 
170 
171 void Foam::vtkPVFoam::updateInfoLagrangian
172 (
173  vtkDataArraySelection* arraySelection
174 )
175 {
176  if (debug)
177  {
178  Info<< "<beg> Foam::vtkPVFoam::updateInfoLagrangian" << nl
179  << " " << dbPtr_->timePath()/cloud::prefix << endl;
180  }
181 
182 
183  // use the db directly since this might be called without a mesh,
184  // but the region must get added back in
185  fileName lagrangianPrefix(cloud::prefix);
186  if (meshRegion_ != polyMesh::defaultRegion)
187  {
188  lagrangianPrefix = meshRegion_/cloud::prefix;
189  }
190 
191  arrayRangeLagrangian_.reset(arraySelection->GetNumberOfArrays());
192 
193  // Generate a list of lagrangian clouds across all times
194  HashSet<fileName> cloudDirs;
195  instantList times = dbPtr_().times();
196  forAll(times, timei)
197  {
198  cloudDirs +=
199  readDir
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  Info<< "<end> Foam::vtkPVFoam::updateInfoLagrangian" << endl;
222  }
223 }
224 
225 
226 void Foam::vtkPVFoam::updateInfoPatches
227 (
228  vtkDataArraySelection* arraySelection,
229  stringList& enabledEntries
230 )
231 {
232  if (debug)
233  {
234  Info<< "<beg> Foam::vtkPVFoam::updateInfoPatches"
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 
253  for
254  (
255  HashTable<labelList, word>::const_iterator iter = groups.begin();
256  iter != groups.end();
257  ++iter
258  )
259  {
260  const word& groupName = iter.key();
261  const labelList& patchIDs = iter();
262 
263  label nFaces = 0;
264  forAll(patchIDs, i)
265  {
266  nFaces += patches[patchIDs[i]].size();
267  }
268 
269  // Valid patch if nFace > 0 - add patch to GUI list
270  if (nFaces)
271  {
272  string vtkGrpName = groupName + " - group";
273  arraySelection->AddArray(vtkGrpName.c_str());
274 
275  ++nPatches;
276 
277  if (enabledEntriesSet.found(vtkGrpName))
278  {
279  if (!reader_->GetShowGroupsOnly())
280  {
281  enabledEntriesSet.erase(vtkGrpName);
282  forAll(patchIDs, i)
283  {
284  const polyPatch& pp = patches[patchIDs[i]];
285  if (pp.size())
286  {
287  string vtkPatchName = pp.name() + " - patch";
288  enabledEntriesSet.insert(vtkPatchName);
289  }
290  }
291  }
292  }
293  }
294  }
295 
296 
297  // Add patches
298  // ~~~~~~~~~~~
299 
300  if (!reader_->GetShowGroupsOnly())
301  {
302  forAll(patches, patchi)
303  {
304  const polyPatch& pp = patches[patchi];
305 
306  if (pp.size())
307  {
308  // Add patch to GUI list
309  arraySelection->AddArray
310  (
311  (pp.name() + " - patch").c_str()
312  );
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  IOobject 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.typeHeaderOk<polyBoundaryMesh>(true))
341  {
342  polyBoundaryMeshEntries patchEntries(ioObj);
343 
344 
345  // Read patches and determine sizes
346  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
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] = readLabel(patchDict.lookup("nFaces"));
356  names[patchi] = patchEntries[patchi].keyword();
357  }
358 
359 
360  // Add (non-zero) patch groups to the list of mesh parts
361  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
362 
363  HashTable<labelList, word> groups(patchEntries.size());
364 
365  forAll(patchEntries, patchi)
366  {
367  const dictionary& patchDict = patchEntries[patchi].dict();
368 
369  wordList groupNames;
370  patchDict.readIfPresent("inGroups", groupNames);
371 
372  forAll(groupNames, groupI)
373  {
374  HashTable<labelList, word>::iterator iter = groups.find
375  (
376  groupNames[groupI]
377  );
378  if (iter != groups.end())
379  {
380  iter().append(patchi);
381  }
382  else
383  {
384  groups.insert(groupNames[groupI], labelList(1, patchi));
385  }
386  }
387  }
388 
389  for
390  (
392  groups.begin();
393  iter != groups.end();
394  ++iter
395  )
396  {
397  const word& groupName = iter.key();
398  const labelList& patchIDs = iter();
399 
400  label nFaces = 0;
401  forAll(patchIDs, i)
402  {
403  nFaces += sizes[patchIDs[i]];
404  }
405 
406  // Valid patch if nFace > 0 - add patch to GUI list
407  if (nFaces)
408  {
409  string vtkGrpName = groupName + " - group";
410  arraySelection->AddArray(vtkGrpName.c_str());
411 
412  ++nPatches;
413 
414  if (enabledEntriesSet.found(vtkGrpName))
415  {
416  if (!reader_->GetShowGroupsOnly())
417  {
418  enabledEntriesSet.erase(vtkGrpName);
419  forAll(patchIDs, i)
420  {
421  if (sizes[patchIDs[i]])
422  {
423  string vtkPatchName =
424  names[patchIDs[i]] + " - patch";
425  enabledEntriesSet.insert(vtkPatchName);
426  }
427  }
428  }
429  }
430  }
431  }
432 
433 
434  // Add (non-zero) patches to the list of mesh parts
435  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
436 
437  if (!reader_->GetShowGroupsOnly())
438  {
439  forAll(names, patchi)
440  {
441  // Valid patch if nFace > 0 - add patch to GUI list
442  if (sizes[patchi])
443  {
444  arraySelection->AddArray
445  (
446  (names[patchi] + " - patch").c_str()
447  );
448 
449  ++nPatches;
450  }
451  }
452  }
453  }
454  }
455  arrayRangePatches_ += nPatches;
456 
457  // Update enabled entries in case of group selection
458  enabledEntries = enabledEntriesSet.toc();
459 
460  if (debug)
461  {
462  // just for debug info
463  getSelectedArrayEntries(arraySelection);
464 
465  Info<< "<end> Foam::vtkPVFoam::updateInfoPatches" << endl;
466  }
467 }
468 
469 
470 void Foam::vtkPVFoam::updateInfoZones
471 (
472  vtkDataArraySelection* arraySelection
473 )
474 {
475  if (!reader_->GetIncludeZones())
476  {
477  return;
478  }
479 
480  if (debug)
481  {
482  Info<< "<beg> Foam::vtkPVFoam::updateInfoZones"
483  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "]" << endl;
484  }
485 
486  wordList namesLst;
487 
488  //
489  // cellZones information
490  // ~~~~~~~~~~~~~~~~~~~~~
491  if (meshPtr_)
492  {
493  namesLst = getZoneNames(meshPtr_->cellZones());
494  }
495  else
496  {
497  namesLst = getZoneNames("cellZones");
498  }
499 
500  arrayRangeCellZones_.reset(arraySelection->GetNumberOfArrays());
501  forAll(namesLst, elemI)
502  {
503  arraySelection->AddArray
504  (
505  (namesLst[elemI] + " - cellZone").c_str()
506  );
507  }
508  arrayRangeCellZones_ += namesLst.size();
509 
510 
511  //
512  // faceZones information
513  // ~~~~~~~~~~~~~~~~~~~~~
514  if (meshPtr_)
515  {
516  namesLst = getZoneNames(meshPtr_->faceZones());
517  }
518  else
519  {
520  namesLst = getZoneNames("faceZones");
521  }
522 
523  arrayRangeFaceZones_.reset(arraySelection->GetNumberOfArrays());
524  forAll(namesLst, elemI)
525  {
526  arraySelection->AddArray
527  (
528  (namesLst[elemI] + " - faceZone").c_str()
529  );
530  }
531  arrayRangeFaceZones_ += namesLst.size();
532 
533 
534  //
535  // pointZones information
536  // ~~~~~~~~~~~~~~~~~~~~~~
537  if (meshPtr_)
538  {
539  namesLst = getZoneNames(meshPtr_->pointZones());
540  }
541  else
542  {
543  namesLst = getZoneNames("pointZones");
544  }
545 
546  arrayRangePointZones_.reset(arraySelection->GetNumberOfArrays());
547  forAll(namesLst, elemI)
548  {
549  arraySelection->AddArray
550  (
551  (namesLst[elemI] + " - pointZone").c_str()
552  );
553  }
554  arrayRangePointZones_ += namesLst.size();
555 
556  if (debug)
557  {
558  // just for debug info
559  getSelectedArrayEntries(arraySelection);
560 
561  Info<< "<end> Foam::vtkPVFoam::updateInfoZones" << endl;
562  }
563 }
564 
565 
566 void Foam::vtkPVFoam::updateInfoSets
567 (
568  vtkDataArraySelection* arraySelection
569 )
570 {
571  if (!reader_->GetIncludeSets())
572  {
573  return;
574  }
575 
576  if (debug)
577  {
578  Info<< "<beg> Foam::vtkPVFoam::updateInfoSets" << endl;
579  }
580 
581  // Add names of sets. Search for last time directory with a sets
582  // subdirectory. Take care not to search beyond the last mesh.
583 
584  word facesInstance = dbPtr_().findInstance
585  (
586  meshDir_,
587  "faces",
589  );
590 
591  word setsInstance = dbPtr_().findInstance
592  (
593  meshDir_/"sets",
594  word::null,
596  facesInstance
597  );
598 
599  IOobjectList objects(dbPtr_(), setsInstance, meshDir_/"sets");
600 
601  if (debug)
602  {
603  Info<< " Foam::vtkPVFoam::updateInfoSets read "
604  << objects.names() << " from " << setsInstance << endl;
605  }
606 
607 
608  arrayRangeCellSets_.reset(arraySelection->GetNumberOfArrays());
609  arrayRangeCellSets_ += addToSelection<cellSet>
610  (
611  arraySelection,
612  objects,
613  " - cellSet"
614  );
615 
616  arrayRangeFaceSets_.reset(arraySelection->GetNumberOfArrays());
617  arrayRangeFaceSets_ += addToSelection<faceSet>
618  (
619  arraySelection,
620  objects,
621  " - faceSet"
622  );
623 
624  arrayRangePointSets_.reset(arraySelection->GetNumberOfArrays());
625  arrayRangePointSets_ += addToSelection<pointSet>
626  (
627  arraySelection,
628  objects,
629  " - pointSet"
630  );
631 
632  if (debug)
633  {
634  // just for debug info
635  getSelectedArrayEntries(arraySelection);
636 
637  Info<< "<end> Foam::vtkPVFoam::updateInfoSets" << endl;
638  }
639 }
640 
641 
642 void Foam::vtkPVFoam::updateInfoLagrangianFields()
643 {
644  if (debug)
645  {
646  Info<< "<beg> Foam::vtkPVFoam::updateInfoLagrangianFields"
647  << endl;
648  }
649 
650  vtkDataArraySelection* fieldSelection =
651  reader_->GetLagrangianFieldSelection();
652 
653  // preserve the enabled selections
654  stringList enabledEntries = getSelectedArrayEntries(fieldSelection);
655  fieldSelection->RemoveAllArrays();
656 
657  // use the db directly since this might be called without a mesh,
658  // but the region must get added back in
659  fileName lagrangianPrefix(cloud::prefix);
660  if (meshRegion_ != polyMesh::defaultRegion)
661  {
662  lagrangianPrefix = meshRegion_/cloud::prefix;
663  }
664 
665  // Add the available fields from all clouds and all time directories.
666  // Differing sets of fields from multiple clouds get combined into a single
667  // set. ParaView will display "(partial)" after field names that only apply
668  // to some of the clouds.
669  const arrayRange& range = arrayRangeLagrangian_;
670  for (label partId = range.start(); partId < range.end(); ++ partId)
671  {
672  const instantList times = dbPtr_().times();
673  forAll(times, timei)
674  {
675  IOobjectList objects
676  (
677  dbPtr_(),
678  times[timei].name(),
679  lagrangianPrefix/getPartName(partId)
680  );
681 
682  addToSelection<IOField<label>>(fieldSelection, objects);
683  addToSelection<IOField<scalar>>(fieldSelection, objects);
684  addToSelection<IOField<vector>>(fieldSelection, objects);
685  addToSelection<IOField<sphericalTensor>>(fieldSelection, objects);
686  addToSelection<IOField<symmTensor>>(fieldSelection, objects);
687  addToSelection<IOField<tensor>>(fieldSelection, objects);
688  }
689  }
690 
691  // restore the enabled selections
692  setSelectedArrayEntries(fieldSelection, enabledEntries);
693 
694  if (debug)
695  {
696  Info<< "<end> Foam::vtkPVFoam::updateInfoLagrangianFields" << endl;
697  }
698 }
699 
700 
701 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:402
List< instant > instantList
List of instants.
Definition: instantList.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:253
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
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
static const char nl
Definition: Ostream.H:262
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
const bool writeData(readBool(pdfDictionary.lookup("writeData")))
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
List< word > wordList
A List of words.
Definition: fileName.H:54
List< string > stringList
A List of strings.
Definition: stringList.H:50
label patchi
fileNameList readDir(const fileName &, const fileName::Type=fileName::FILE, const bool filtergz=true, const bool followLink=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:641
messageStream Info
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.