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-2016 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("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.headerOk())
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  // Search for list of lagrangian objects for this time
192  fileNameList cloudDirs
193  (
194  readDir(dbPtr_->timePath()/lagrangianPrefix, fileName::DIRECTORY)
195  );
196 
197  arrayRangeLagrangian_.reset(arraySelection->GetNumberOfArrays());
198 
199  int nClouds = 0;
200  forAll(cloudDirs, cloudI)
201  {
202  // Add cloud to GUI list
203  arraySelection->AddArray
204  (
205  (cloudDirs[cloudI] + " - lagrangian").c_str()
206  );
207 
208  ++nClouds;
209  }
210  arrayRangeLagrangian_ += nClouds;
211 
212  if (debug)
213  {
214  // just for debug info
215  getSelectedArrayEntries(arraySelection);
216 
217  Info<< "<end> Foam::vtkPVFoam::updateInfoLagrangian" << endl;
218  }
219 }
220 
221 
222 void Foam::vtkPVFoam::updateInfoPatches
223 (
224  vtkDataArraySelection* arraySelection,
225  stringList& enabledEntries
226 )
227 {
228  if (debug)
229  {
230  Info<< "<beg> Foam::vtkPVFoam::updateInfoPatches"
231  << " [meshPtr=" << (meshPtr_ ? "set" : "NULL") << "]" << endl;
232  }
233 
234 
235  HashSet<string> enabledEntriesSet(enabledEntries);
236 
237  arrayRangePatches_.reset(arraySelection->GetNumberOfArrays());
238 
239  int nPatches = 0;
240  if (meshPtr_)
241  {
242  const polyBoundaryMesh& patches = meshPtr_->boundaryMesh();
243  const HashTable<labelList, word>& groups = patches.groupPatchIDs();
244  const wordList allPatchNames = patches.names();
245 
246  // Add patch groups
247  // ~~~~~~~~~~~~~~~~
248 
249  for
250  (
251  HashTable<labelList, word>::const_iterator iter = groups.begin();
252  iter != groups.end();
253  ++iter
254  )
255  {
256  const word& groupName = iter.key();
257  const labelList& patchIDs = iter();
258 
259  label nFaces = 0;
260  forAll(patchIDs, i)
261  {
262  nFaces += patches[patchIDs[i]].size();
263  }
264 
265  // Valid patch if nFace > 0 - add patch to GUI list
266  if (nFaces)
267  {
268  string vtkGrpName = groupName + " - group";
269  arraySelection->AddArray(vtkGrpName.c_str());
270 
271  ++nPatches;
272 
273  if (enabledEntriesSet.found(vtkGrpName))
274  {
275  if (!reader_->GetShowGroupsOnly())
276  {
277  enabledEntriesSet.erase(vtkGrpName);
278  forAll(patchIDs, i)
279  {
280  const polyPatch& pp = patches[patchIDs[i]];
281  if (pp.size())
282  {
283  string vtkPatchName = pp.name() + " - patch";
284  enabledEntriesSet.insert(vtkPatchName);
285  }
286  }
287  }
288  }
289  }
290  }
291 
292 
293  // Add patches
294  // ~~~~~~~~~~~
295 
296  if (!reader_->GetShowGroupsOnly())
297  {
298  forAll(patches, patchi)
299  {
300  const polyPatch& pp = patches[patchi];
301 
302  if (pp.size())
303  {
304  // Add patch to GUI list
305  arraySelection->AddArray
306  (
307  (pp.name() + " - patch").c_str()
308  );
309 
310  ++nPatches;
311  }
312  }
313  }
314  }
315  else
316  {
317  // mesh not loaded - read from file
318  // but this could fail if we've supplied a bad region name
319  IOobject ioObj
320  (
321  "boundary",
322  dbPtr_().findInstance
323  (
324  meshDir_,
325  "boundary",
327  ),
328  meshDir_,
329  dbPtr_(),
332  false
333  );
334 
335  // this should only ever fail if the mesh region doesn't exist
336  if (ioObj.headerOk())
337  {
338  polyBoundaryMeshEntries patchEntries(ioObj);
339 
340 
341  // Read patches and determine sizes
342  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
343 
344  wordList names(patchEntries.size());
345  labelList sizes(patchEntries.size());
346 
347  forAll(patchEntries, patchi)
348  {
349  const dictionary& patchDict = patchEntries[patchi].dict();
350 
351  sizes[patchi] = readLabel(patchDict.lookup("nFaces"));
352  names[patchi] = patchEntries[patchi].keyword();
353  }
354 
355 
356  // Add (non-zero) patch groups to the list of mesh parts
357  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
358 
359  HashTable<labelList, word> groups(patchEntries.size());
360 
361  forAll(patchEntries, patchi)
362  {
363  const dictionary& patchDict = patchEntries[patchi].dict();
364 
365  wordList groupNames;
366  patchDict.readIfPresent("inGroups", groupNames);
367 
368  forAll(groupNames, groupI)
369  {
370  HashTable<labelList, word>::iterator iter = groups.find
371  (
372  groupNames[groupI]
373  );
374  if (iter != groups.end())
375  {
376  iter().append(patchi);
377  }
378  else
379  {
380  groups.insert(groupNames[groupI], labelList(1, patchi));
381  }
382  }
383  }
384 
385  for
386  (
388  groups.begin();
389  iter != groups.end();
390  ++iter
391  )
392  {
393  const word& groupName = iter.key();
394  const labelList& patchIDs = iter();
395 
396  label nFaces = 0;
397  forAll(patchIDs, i)
398  {
399  nFaces += sizes[patchIDs[i]];
400  }
401 
402  // Valid patch if nFace > 0 - add patch to GUI list
403  if (nFaces)
404  {
405  string vtkGrpName = groupName + " - group";
406  arraySelection->AddArray(vtkGrpName.c_str());
407 
408  ++nPatches;
409 
410  if (enabledEntriesSet.found(vtkGrpName))
411  {
412  if (!reader_->GetShowGroupsOnly())
413  {
414  enabledEntriesSet.erase(vtkGrpName);
415  forAll(patchIDs, i)
416  {
417  if (sizes[patchIDs[i]])
418  {
419  string vtkPatchName =
420  names[patchIDs[i]] + " - patch";
421  enabledEntriesSet.insert(vtkPatchName);
422  }
423  }
424  }
425  }
426  }
427  }
428 
429 
430  // Add (non-zero) patches to the list of mesh parts
431  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
432 
433  if (!reader_->GetShowGroupsOnly())
434  {
435  forAll(names, patchi)
436  {
437  // Valid patch if nFace > 0 - add patch to GUI list
438  if (sizes[patchi])
439  {
440  arraySelection->AddArray
441  (
442  (names[patchi] + " - patch").c_str()
443  );
444 
445  ++nPatches;
446  }
447  }
448  }
449  }
450  }
451  arrayRangePatches_ += nPatches;
452 
453  // Update enabled entries in case of group selection
454  enabledEntries = enabledEntriesSet.toc();
455 
456  if (debug)
457  {
458  // just for debug info
459  getSelectedArrayEntries(arraySelection);
460 
461  Info<< "<end> Foam::vtkPVFoam::updateInfoPatches" << endl;
462  }
463 }
464 
465 
466 void Foam::vtkPVFoam::updateInfoZones
467 (
468  vtkDataArraySelection* arraySelection
469 )
470 {
471  if (!reader_->GetIncludeZones())
472  {
473  return;
474  }
475 
476  if (debug)
477  {
478  Info<< "<beg> Foam::vtkPVFoam::updateInfoZones"
479  << " [meshPtr=" << (meshPtr_ ? "set" : "NULL") << "]" << endl;
480  }
481 
482  wordList namesLst;
483 
484  //
485  // cellZones information
486  // ~~~~~~~~~~~~~~~~~~~~~
487  if (meshPtr_)
488  {
489  namesLst = getZoneNames(meshPtr_->cellZones());
490  }
491  else
492  {
493  namesLst = getZoneNames("cellZones");
494  }
495 
496  arrayRangeCellZones_.reset(arraySelection->GetNumberOfArrays());
497  forAll(namesLst, elemI)
498  {
499  arraySelection->AddArray
500  (
501  (namesLst[elemI] + " - cellZone").c_str()
502  );
503  }
504  arrayRangeCellZones_ += namesLst.size();
505 
506 
507  //
508  // faceZones information
509  // ~~~~~~~~~~~~~~~~~~~~~
510  if (meshPtr_)
511  {
512  namesLst = getZoneNames(meshPtr_->faceZones());
513  }
514  else
515  {
516  namesLst = getZoneNames("faceZones");
517  }
518 
519  arrayRangeFaceZones_.reset(arraySelection->GetNumberOfArrays());
520  forAll(namesLst, elemI)
521  {
522  arraySelection->AddArray
523  (
524  (namesLst[elemI] + " - faceZone").c_str()
525  );
526  }
527  arrayRangeFaceZones_ += namesLst.size();
528 
529 
530  //
531  // pointZones information
532  // ~~~~~~~~~~~~~~~~~~~~~~
533  if (meshPtr_)
534  {
535  namesLst = getZoneNames(meshPtr_->pointZones());
536  }
537  else
538  {
539  namesLst = getZoneNames("pointZones");
540  }
541 
542  arrayRangePointZones_.reset(arraySelection->GetNumberOfArrays());
543  forAll(namesLst, elemI)
544  {
545  arraySelection->AddArray
546  (
547  (namesLst[elemI] + " - pointZone").c_str()
548  );
549  }
550  arrayRangePointZones_ += namesLst.size();
551 
552  if (debug)
553  {
554  // just for debug info
555  getSelectedArrayEntries(arraySelection);
556 
557  Info<< "<end> Foam::vtkPVFoam::updateInfoZones" << endl;
558  }
559 }
560 
561 
562 void Foam::vtkPVFoam::updateInfoSets
563 (
564  vtkDataArraySelection* arraySelection
565 )
566 {
567  if (!reader_->GetIncludeSets())
568  {
569  return;
570  }
571 
572  if (debug)
573  {
574  Info<< "<beg> Foam::vtkPVFoam::updateInfoSets" << endl;
575  }
576 
577  // Add names of sets. Search for last time directory with a sets
578  // subdirectory. Take care not to search beyond the last mesh.
579 
580  word facesInstance = dbPtr_().findInstance
581  (
582  meshDir_,
583  "faces",
585  );
586 
587  word setsInstance = dbPtr_().findInstance
588  (
589  meshDir_/"sets",
590  word::null,
592  facesInstance
593  );
594 
595  IOobjectList objects(dbPtr_(), setsInstance, meshDir_/"sets");
596 
597  if (debug)
598  {
599  Info<< " Foam::vtkPVFoam::updateInfoSets read "
600  << objects.names() << " from " << setsInstance << endl;
601  }
602 
603 
604  arrayRangeCellSets_.reset(arraySelection->GetNumberOfArrays());
605  arrayRangeCellSets_ += addToSelection<cellSet>
606  (
607  arraySelection,
608  objects,
609  " - cellSet"
610  );
611 
612  arrayRangeFaceSets_.reset(arraySelection->GetNumberOfArrays());
613  arrayRangeFaceSets_ += addToSelection<faceSet>
614  (
615  arraySelection,
616  objects,
617  " - faceSet"
618  );
619 
620  arrayRangePointSets_.reset(arraySelection->GetNumberOfArrays());
621  arrayRangePointSets_ += addToSelection<pointSet>
622  (
623  arraySelection,
624  objects,
625  " - pointSet"
626  );
627 
628  if (debug)
629  {
630  // just for debug info
631  getSelectedArrayEntries(arraySelection);
632 
633  Info<< "<end> Foam::vtkPVFoam::updateInfoSets" << endl;
634  }
635 }
636 
637 
638 void Foam::vtkPVFoam::updateInfoLagrangianFields()
639 {
640  if (debug)
641  {
642  Info<< "<beg> Foam::vtkPVFoam::updateInfoLagrangianFields"
643  << endl;
644  }
645 
646  vtkDataArraySelection* fieldSelection =
647  reader_->GetLagrangianFieldSelection();
648 
649  // preserve the enabled selections
650  stringList enabledEntries = getSelectedArrayEntries(fieldSelection);
651  fieldSelection->RemoveAllArrays();
652 
653  // TODO - currently only get fields from ONE cloud
654  // have to decide if the second set of fields get mixed in
655  // or dealt with separately
656 
657  const arrayRange& range = arrayRangeLagrangian_;
658  if (range.empty())
659  {
660  return;
661  }
662 
663  int partId = range.start();
664  word cloudName = getPartName(partId);
665 
666  // use the db directly since this might be called without a mesh,
667  // but the region must get added back in
668  fileName lagrangianPrefix(cloud::prefix);
669  if (meshRegion_ != polyMesh::defaultRegion)
670  {
671  lagrangianPrefix = meshRegion_/cloud::prefix;
672  }
673 
674  IOobjectList objects
675  (
676  dbPtr_(),
677  dbPtr_().timeName(),
678  lagrangianPrefix/cloudName
679  );
680 
681  addToSelection<IOField<label>>
682  (
683  fieldSelection,
684  objects
685  );
686  addToSelection<IOField<scalar>>
687  (
688  fieldSelection,
689  objects
690  );
691  addToSelection<IOField<vector>>
692  (
693  fieldSelection,
694  objects
695  );
696  addToSelection<IOField<sphericalTensor>>
697  (
698  fieldSelection,
699 
700  objects
701  );
702  addToSelection<IOField<symmTensor>>
703  (
704  fieldSelection,
705  objects
706  );
707  addToSelection<IOField<tensor>>
708  (
709  fieldSelection,
710  objects
711  );
712 
713  // restore the enabled selections
714  setSelectedArrayEntries(fieldSelection, enabledEntries);
715 
716  if (debug)
717  {
718  Info<< "<end> Foam::vtkPVFoam::updateInfoLagrangianFields - "
719  << "lagrangian objects.size() = " << objects.size() << endl;
720  }
721 }
722 
723 
724 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:402
#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:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
stringList getSelectedArrayEntries(vtkDataArraySelection *)
Retrieve the current selections.
friend class const_iterator
Declare friendship with the const_iterator.
Definition: HashTable.H:192
static const word null
An empty word.
Definition: word.H:77
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
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")))
List< word > wordList
A List of words.
Definition: fileName.H:54
List< string > stringList
A List of strings.
Definition: stringList.H:50
label patchi
void setSelectedArrayEntries(vtkDataArraySelection *, const stringList &)
Set selection(s)
messageStream Info
friend class iterator
Declare friendship with the iterator.
Definition: HashTable.H:189
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
fileNameList readDir(const fileName &, const fileName::Type=fileName::FILE, const bool filtergz=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:527
Namespace for OpenFOAM.