vtkPVFoam.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-2018 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 #include "vtkPVFoamReader.h"
28 
29 // OpenFOAM includes
30 #include "fvMesh.H"
31 #include "Time.H"
32 #include "patchZones.H"
33 #include "collatedFileOperation.H"
34 #include "etcFiles.H"
35 
36 // VTK includes
37 #include "vtkDataArraySelection.h"
38 #include "vtkMultiBlockDataSet.h"
39 #include "vtkRenderer.h"
40 #include "vtkTextActor.h"
41 #include "vtkTextProperty.h"
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47  defineTypeNameAndDebug(vtkPVFoam, 0);
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
55 
56 void Foam::vtkPVFoam::resetCounters()
57 {
58  // Reset array range information (ids and sizes)
59  arrayRangeVolume_.reset();
60  arrayRangePatches_.reset();
61  arrayRangeLagrangian_.reset();
62  arrayRangeCellZones_.reset();
63  arrayRangeFaceZones_.reset();
64  arrayRangePointZones_.reset();
65  arrayRangeCellSets_.reset();
66  arrayRangeFaceSets_.reset();
67  arrayRangePointSets_.reset();
68 }
69 
70 
71 void Foam::vtkPVFoam::reduceMemory()
72 {
73  forAll(regionPolyDecomp_, i)
74  {
75  regionPolyDecomp_[i].clear();
76  }
77 
78  forAll(zonePolyDecomp_, i)
79  {
80  zonePolyDecomp_[i].clear();
81  }
82 
83  forAll(csetPolyDecomp_, i)
84  {
85  csetPolyDecomp_[i].clear();
86  }
87 
88  if (!reader_->GetCacheMesh())
89  {
90  delete meshPtr_;
91  meshPtr_ = nullptr;
92  }
93 }
94 
95 
96 int Foam::vtkPVFoam::setTime(int nRequest, const double requestTimes[])
97 {
98  Time& runTime = dbPtr_();
99 
100  // Get times list. Flush first to force refresh.
101  fileHandler().flush();
102  instantList Times = runTime.times();
103 
104  int nearestIndex = timeIndex_;
105  for (int requestI = 0; requestI < nRequest; ++requestI)
106  {
107  int index = Time::findClosestTimeIndex(Times, requestTimes[requestI]);
108  if (index >= 0 && index != timeIndex_)
109  {
110  nearestIndex = index;
111  break;
112  }
113  }
114 
115  if (nearestIndex < 0)
116  {
117  nearestIndex = 0;
118  }
119 
120  if (debug)
121  {
122  Info<< "<beg> Foam::vtkPVFoam::setTime(";
123  for (int requestI = 0; requestI < nRequest; ++requestI)
124  {
125  if (requestI)
126  {
127  Info<< ", ";
128  }
129 
130  Info<< requestTimes[requestI];
131  }
132  Info<< ") - previousIndex = " << timeIndex_
133  << ", nearestIndex = " << nearestIndex << endl;
134  }
135 
136 
137  // see what has changed
138  if (timeIndex_ != nearestIndex)
139  {
140  timeIndex_ = nearestIndex;
141  runTime.setTime(Times[nearestIndex], nearestIndex);
142 
143  // the fields change each time
144  fieldsChanged_ = true;
145 
146  if (meshPtr_)
147  {
148  if (meshPtr_->readUpdate() != polyMesh::UNCHANGED)
149  {
150  meshChanged_ = true;
151  }
152  }
153  else
154  {
155  meshChanged_ = true;
156  }
157  }
158 
159  if (debug)
160  {
161  Info<< "<end> Foam::vtkPVFoam::setTime() - selectedTime="
162  << Times[nearestIndex].name() << " index=" << timeIndex_
163  << "/" << Times.size()
164  << " meshChanged=" << Switch(meshChanged_)
165  << " fieldsChanged=" << Switch(fieldsChanged_) << endl;
166  }
167 
168  return nearestIndex;
169 }
170 
171 
172 void Foam::vtkPVFoam::updateMeshPartsStatus()
173 {
174  if (debug)
175  {
176  Info<< "<beg> Foam::vtkPVFoam::updateMeshPartsStatus" << endl;
177  }
178 
179  vtkDataArraySelection* selection = reader_->GetPartSelection();
180  label nElem = selection->GetNumberOfArrays();
181 
182  if (partStatus_.size() != nElem)
183  {
184  partStatus_.setSize(nElem);
185  partStatus_ = false;
186  meshChanged_ = true;
187  }
188 
189  // this needs fixing if we wish to re-use the datasets
190  partDataset_.setSize(nElem);
191  partDataset_ = -1;
192 
193  // Read the selected mesh parts (zones, patches ...) and add to list
194  forAll(partStatus_, partId)
195  {
196  const int setting = selection->GetArraySetting(partId);
197 
198  if (partStatus_[partId] != setting)
199  {
200  partStatus_[partId] = setting;
201  meshChanged_ = true;
202  }
203 
204  if (debug)
205  {
206  Info<< " part[" << partId << "] = "
207  << partStatus_[partId]
208  << " : " << selection->GetArrayName(partId) << endl;
209  }
210  }
211  if (debug)
212  {
213  Info<< "<end> Foam::vtkPVFoam::updateMeshPartsStatus" << endl;
214  }
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
219 
221 (
222  const char* const vtkFileName,
223  vtkPVFoamReader* reader
224 )
225 :
226  reader_(reader),
227  dbPtr_(nullptr),
228  meshPtr_(nullptr),
229  meshRegion_(polyMesh::defaultRegion),
230  meshDir_(polyMesh::meshSubDir),
231  timeIndex_(-1),
232  meshChanged_(true),
233  fieldsChanged_(true),
234  arrayRangeVolume_("unzoned"),
235  arrayRangePatches_("patches"),
236  arrayRangeLagrangian_("lagrangian"),
237  arrayRangeCellZones_("cellZone"),
238  arrayRangeFaceZones_("faceZone"),
239  arrayRangePointZones_("pointZone"),
240  arrayRangeCellSets_("cellSet"),
241  arrayRangeFaceSets_("faceSet"),
242  arrayRangePointSets_("pointSet")
243 {
244  if (debug)
245  {
246  Info<< "Foam::vtkPVFoam::vtkPVFoam - " << vtkFileName << endl;
247  printMemory();
248  }
249 
250  fileName FileName(vtkFileName);
251 
252  // avoid argList and get rootPath/caseName directly from the file
253  fileName fullCasePath(FileName.path());
254 
255  if (!isDir(fullCasePath))
256  {
257  return;
258  }
259  if (fullCasePath == ".")
260  {
261  fullCasePath = cwd();
262  }
263 
264 
265  if (fullCasePath.name().find("processors", 0) == 0)
266  {
267  // FileName e.g. "cavity/processors256/processor1.OpenFOAM
268  // Remove the processors section so it goes into processorDDD
269  // checking below.
270  fullCasePath = fullCasePath.path()/fileName(FileName.name()).lessExt();
271  }
272 
273 
274  if (fullCasePath.name().find("processor", 0) == 0)
275  {
276  // Give filehandler opportunity to analyse number of processors
277  (void)fileHandler().filePath(fullCasePath);
278 
279  const fileName globalCase = fullCasePath.path();
280 
281  setEnv("FOAM_CASE", globalCase, true);
282  setEnv("FOAM_CASENAME", globalCase.name(), true);
283  }
284  else
285  {
286  setEnv("FOAM_CASE", fullCasePath, true);
287  setEnv("FOAM_CASENAME", fullCasePath.name(), true);
288  }
289 
290  // look for 'case{region}.OpenFOAM'
291  // could be stringent and insist the prefix match the directory name...
292  // Note: cannot use fileName::name() due to the embedded '{}'
293  string caseName(FileName.lessExt());
294  string::size_type beg = caseName.find_last_of("/{");
295  string::size_type end = caseName.find('}', beg);
296 
297  if
298  (
299  beg != string::npos && caseName[beg] == '{'
300  && end != string::npos && end == caseName.size()-1
301  )
302  {
303  meshRegion_ = caseName.substr(beg+1, end-beg-1);
304 
305  // some safety
306  if (meshRegion_.empty())
307  {
308  meshRegion_ = polyMesh::defaultRegion;
309  }
310 
311  if (meshRegion_ != polyMesh::defaultRegion)
312  {
313  meshDir_ = meshRegion_/polyMesh::meshSubDir;
314  }
315  }
316 
317  if (debug)
318  {
319  Info<< "fullCasePath=" << fullCasePath << nl
320  << "FOAM_CASE=" << getEnv("FOAM_CASE") << nl
321  << "FOAM_CASENAME=" << getEnv("FOAM_CASENAME") << nl
322  << "region=" << meshRegion_ << endl;
323  }
324 
325  // Pre-load any libraries
326  dlLibraryTable dlTable;
327 
328  string libsString(getEnv("FOAM_LIBS"));
329  if (!libsString.empty())
330  {
331  IStringStream is(libsString);
332  fileNameList libNames(is);
333  forAll(libNames, i)
334  {
335  dlTable.open(libNames[i]);
336  }
337  }
338 
339  // Create time object
340  dbPtr_.reset
341  (
342  new Time
343  (
344  Time::controlDictName,
345  fileName(fullCasePath.path()),
346  fileName(fullCasePath.name())
347  )
348  );
349 
350  dbPtr_().functionObjects().off();
351 
352  fileNameList configDictFiles = findEtcFiles("paraFoam", false);
353  forAllReverse(configDictFiles, cdfi)
354  {
355  configDict_.merge(dictionary(IFstream(configDictFiles[cdfi])()));
356  }
357 
358  updateInfo();
359 }
360 
361 
362 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
363 
365 {
366  if (debug)
367  {
368  Info<< "<end> Foam::vtkPVFoam::~vtkPVFoam" << endl;
369  }
370 
371  delete meshPtr_;
372 }
373 
374 
375 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
376 
378 {
379  if (debug)
380  {
381  Info<< "<beg> Foam::vtkPVFoam::updateInfo"
382  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "] timeIndex="
383  << timeIndex_ << endl;
384  }
385 
386  resetCounters();
387 
388  vtkDataArraySelection* partSelection = reader_->GetPartSelection();
389 
390  // there are two ways to ensure we have the correct list of parts:
391  // 1. remove everything and then set particular entries 'on'
392  // 2. build a 'char **' list and call SetArraysWithDefault()
393  //
394  // Nr. 2 has the potential advantage of not touching the modification
395  // time of the vtkDataArraySelection, but the qt/paraview proxy
396  // layer doesn't care about that anyhow.
397 
398  // enable 'internalMesh' on the first call
399  // or preserve the enabled selections
400  stringList enabledEntries;
401  bool first = !partSelection->GetNumberOfArrays() && !meshPtr_;
402  if (first)
403  {
404  enabledEntries.setSize(1);
405  enabledEntries[0] = "internalMesh";
406  }
407  else
408  {
409  enabledEntries = getSelectedArrayEntries(partSelection);
410  }
411 
412  // Clear current mesh parts list
413  partSelection->RemoveAllArrays();
414 
415  // Update mesh parts list - add Lagrangian at the bottom
416  updateInfoInternalMesh(partSelection);
417  updateInfoPatches(partSelection, enabledEntries, first);
418  updateInfoSets(partSelection);
419  updateInfoZones(partSelection);
420  updateInfoLagrangian(partSelection);
421 
422  // restore the enabled selections
423  setSelectedArrayEntries(partSelection, enabledEntries);
424 
425  if (meshChanged_)
426  {
427  fieldsChanged_ = true;
428  }
429 
430  // Update volume, point and lagrangian fields
431  updateInfoFields<fvPatchField, volMesh>
432  (
433  reader_->GetVolFieldSelection()
434  );
435  updateInfoFields<pointPatchField, pointMesh>
436  (
437  reader_->GetPointFieldSelection()
438  );
439  updateInfoLagrangianFields();
440 
441  if (debug)
442  {
443  // just for debug info
444  getSelectedArrayEntries(partSelection);
445  Info<< "<end> Foam::vtkPVFoam::updateInfo" << endl;
446  }
447 
448 }
449 
450 
451 void Foam::vtkPVFoam::updateFoamMesh()
452 {
453  if (debug)
454  {
455  Info<< "<beg> Foam::vtkPVFoam::updateFoamMesh" << endl;
456  printMemory();
457  }
458 
459  if (!reader_->GetCacheMesh())
460  {
461  delete meshPtr_;
462  meshPtr_ = nullptr;
463  }
464 
465  // Check to see if the OpenFOAM mesh has been created
466  if (!meshPtr_)
467  {
468  if (debug)
469  {
470  Info<< "Creating OpenFOAM mesh for region " << meshRegion_
471  << " at time=" << dbPtr_().timeName()
472  << endl;
473 
474  }
475 
476  meshPtr_ = new fvMesh
477  (
478  IOobject
479  (
480  meshRegion_,
481  dbPtr_().timeName(),
482  dbPtr_(),
484  )
485  );
486 
487  meshChanged_ = true;
488  }
489  else
490  {
491  if (debug)
492  {
493  Info<< "Using existing OpenFOAM mesh" << endl;
494  }
495  }
496 
497  if (debug)
498  {
499  Info<< "<end> Foam::vtkPVFoam::updateFoamMesh" << endl;
500  printMemory();
501  }
502 }
503 
504 
506 (
507  vtkMultiBlockDataSet* output,
508  vtkMultiBlockDataSet* lagrangianOutput
509 )
510 {
511  if (debug)
512  {
513  cout<< "<beg> Foam::vtkPVFoam::Update - output with "
514  << output->GetNumberOfBlocks() << " and "
515  << lagrangianOutput->GetNumberOfBlocks() << " blocks\n";
516  output->Print(cout);
517  lagrangianOutput->Print(cout);
518  printMemory();
519  }
520  reader_->UpdateProgress(0.1);
521 
522  // Set up mesh parts selection(s)
523  updateMeshPartsStatus();
524 
525  reader_->UpdateProgress(0.15);
526 
527  // Update the OpenFOAM mesh
528  updateFoamMesh();
529  reader_->UpdateProgress(0.4);
530 
531  // Convert meshes - start port0 at block=0
532  int blockNo = 0;
533 
534  convertMeshVolume(output, blockNo);
535  convertMeshPatches(output, blockNo);
536  reader_->UpdateProgress(0.6);
537 
538  if (reader_->GetIncludeZones())
539  {
540  convertMeshCellZones(output, blockNo);
541  convertMeshFaceZones(output, blockNo);
542  convertMeshPointZones(output, blockNo);
543  reader_->UpdateProgress(0.65);
544  }
545 
546  if (reader_->GetIncludeSets())
547  {
548  convertMeshCellSets(output, blockNo);
549  convertMeshFaceSets(output, blockNo);
550  convertMeshPointSets(output, blockNo);
551  reader_->UpdateProgress(0.7);
552  }
553 
554  convertMeshLagrangian(lagrangianOutput, blockNo);
555 
556  reader_->UpdateProgress(0.8);
557 
558  // Update fields
559  convertVolFields(output);
560  convertPointFields(output);
561  convertLagrangianFields(lagrangianOutput);
562  if (debug)
563  {
564  Info<< "done reader part" << endl;
565  }
566  reader_->UpdateProgress(0.95);
567 
568  meshChanged_ = fieldsChanged_ = false;
569 }
570 
571 
573 {
574  // reclaim some memory
575  reduceMemory();
576  reader_->UpdateProgress(1.0);
577 }
578 
579 
580 double* Foam::vtkPVFoam::findTimes(int& nTimeSteps)
581 {
582  int nTimes = 0;
583  double* tsteps = nullptr;
584 
585  if (dbPtr_.valid())
586  {
587  Time& runTime = dbPtr_();
588  // Get times list. Flush first to force refresh.
589  fileHandler().flush();
590  instantList timeLst = runTime.times();
591 
592  // find the first time for which this mesh appears to exist
593  label timeI = 0;
594  for (; timeI < timeLst.size(); ++timeI)
595  {
596  const word& timeName = timeLst[timeI].name();
597 
598  if
599  (
600  IOobject
601  (
602  "points",
603  timeName,
604  meshDir_,
605  runTime
606  ).typeHeaderOk<pointIOField>(true)
607  )
608  {
609  break;
610  }
611  }
612 
613  nTimes = timeLst.size() - timeI;
614 
615  // skip "constant" time whenever possible
616  if (timeI == 0 && nTimes > 1)
617  {
618  if (timeLst[timeI].name() == runTime.constant())
619  {
620  ++timeI;
621  --nTimes;
622  }
623  }
624 
625 
626  // skip "0/" time if requested and possible
627  if (nTimes > 1 && reader_->GetSkipZeroTime())
628  {
629  if (mag(timeLst[timeI].value()) < small)
630  {
631  ++timeI;
632  --nTimes;
633  }
634  }
635 
636  if (nTimes)
637  {
638  tsteps = new double[nTimes];
639  for (label stepI = 0; stepI < nTimes; ++stepI, ++timeI)
640  {
641  tsteps[stepI] = timeLst[timeI].value();
642  }
643  }
644  }
645  else
646  {
647  if (debug)
648  {
649  cout<< "no valid dbPtr:\n";
650  }
651  }
652 
653  // vector length returned via the parameter
654  nTimeSteps = nTimes;
655 
656  return tsteps;
657 }
658 
659 
661 (
662  vtkRenderer* renderer,
663  const bool show
664 )
665 {
666  if (!meshPtr_)
667  {
668  return;
669  }
670 
671  // always remove old actors first
672 
673  forAll(patchTextActorsPtrs_, patchi)
674  {
675  renderer->RemoveViewProp(patchTextActorsPtrs_[patchi]);
676  patchTextActorsPtrs_[patchi]->Delete();
677  }
678  patchTextActorsPtrs_.clear();
679 
680  if (show)
681  {
682  // get the display patches, strip off any suffix
683  wordHashSet selectedPatches = getSelected
684  (
685  reader_->GetPartSelection(),
686  arrayRangePatches_
687  );
688 
689  if (selectedPatches.empty())
690  {
691  return;
692  }
693 
694  const polyBoundaryMesh& pbMesh = meshPtr_->boundaryMesh();
695 
696  // Find the total number of zones
697  // Each zone will take the patch name
698  // Number of zones per patch ... zero zones should be skipped
699  labelList nZones(pbMesh.size(), 0);
700 
701  // Per global zone number the average face centre position
702  List<DynamicList<point>> zoneCentre(pbMesh.size());
703 
704 
705  // Loop through all patches to determine zones, and centre of each zone
706  forAll(pbMesh, patchi)
707  {
708  const polyPatch& pp = pbMesh[patchi];
709 
710  // Only include the patch if it is selected
711  if (!selectedPatches.found(pp.name()))
712  {
713  continue;
714  }
715 
716  const labelListList& edgeFaces = pp.edgeFaces();
717  const vectorField& n = pp.faceNormals();
718 
719  boolList featEdge(pp.nEdges(), false);
720 
721  forAll(edgeFaces, edgeI)
722  {
723  const labelList& eFaces = edgeFaces[edgeI];
724 
725  if (eFaces.size() == 1)
726  {
727  // Note: could also do ones with > 2 faces but this gives
728  // too many zones for baffles
729  featEdge[edgeI] = true;
730  }
731  else if (mag(n[eFaces[0]] & n[eFaces[1]]) < 0.5)
732  {
733  featEdge[edgeI] = true;
734  }
735  }
736 
737  // Do topological analysis of patch, find disconnected regions
738  patchZones pZones(pp, featEdge);
739 
740  nZones[patchi] = pZones.nZones();
741 
742  labelList zoneNFaces(pZones.nZones(), 0);
743 
744  // Create storage for additional zone centres
745  forAll(zoneNFaces, zoneI)
746  {
747  zoneCentre[patchi].append(Zero);
748  }
749 
750  // Do averaging per individual zone
751  forAll(pp, facei)
752  {
753  label zoneI = pZones[facei];
754  zoneCentre[patchi][zoneI] += pp[facei].centre(pp.points());
755  zoneNFaces[zoneI]++;
756  }
757 
758  forAll(zoneCentre[patchi], zoneI)
759  {
760  zoneCentre[patchi][zoneI] /= zoneNFaces[zoneI];
761  }
762  }
763 
764  // Count number of zones we're actually going to display.
765  // This is truncated to a max per patch
766 
767  const label MAXPATCHZONES = 20;
768 
769  label displayZoneI = 0;
770 
771  forAll(pbMesh, patchi)
772  {
773  displayZoneI += min(MAXPATCHZONES, nZones[patchi]);
774  }
775 
776  if (debug)
777  {
778  Info<< "displayed zone centres = " << displayZoneI << nl
779  << "zones per patch = " << nZones << endl;
780  }
781 
782  // Set the size of the patch labels to max number of zones
783  patchTextActorsPtrs_.setSize(displayZoneI);
784 
785  if (debug)
786  {
787  Info<< "constructing patch labels" << endl;
788  }
789 
790  // Actor index
791  displayZoneI = 0;
792 
793  forAll(pbMesh, patchi)
794  {
795  const polyPatch& pp = pbMesh[patchi];
796 
797  label globalZoneI = 0;
798 
799  // Only selected patches will have a non-zero number of zones
800  label nDisplayZones = min(MAXPATCHZONES, nZones[patchi]);
801  label increment = 1;
802  if (nZones[patchi] >= MAXPATCHZONES)
803  {
804  increment = nZones[patchi]/MAXPATCHZONES;
805  }
806 
807  for (label i = 0; i < nDisplayZones; i++)
808  {
809  if (debug)
810  {
811  Info<< "patch name = " << pp.name() << nl
812  << "anchor = " << zoneCentre[patchi][globalZoneI] << nl
813  << "globalZoneI = " << globalZoneI << endl;
814  }
815 
816  vtkTextActor* txt = vtkTextActor::New();
817 
818  txt->SetInput(pp.name().c_str());
819 
820  // Set text properties
821  vtkTextProperty* tprop = txt->GetTextProperty();
822  tprop->SetFontFamilyToArial();
823  tprop->BoldOff();
824  tprop->ShadowOff();
825  tprop->SetLineSpacing(1.0);
826  tprop->SetFontSize(12);
827  tprop->SetColor(1.0, 0.0, 0.0);
828  tprop->SetJustificationToCentered();
829 
830  // Set text to use 3-D world co-ordinates
831  txt->GetPositionCoordinate()->SetCoordinateSystemToWorld();
832 
833  txt->GetPositionCoordinate()->SetValue
834  (
835  zoneCentre[patchi][globalZoneI].x(),
836  zoneCentre[patchi][globalZoneI].y(),
837  zoneCentre[patchi][globalZoneI].z()
838  );
839 
840  // Add text to each renderer
841  renderer->AddViewProp(txt);
842 
843  // Maintain a list of text labels added so that they can be
844  // removed later
845  patchTextActorsPtrs_[displayZoneI] = txt;
846 
847  globalZoneI += increment;
848  displayZoneI++;
849  }
850  }
851 
852  // Resize the patch names list to the actual number of patch names added
853  patchTextActorsPtrs_.setSize(displayZoneI);
854  }
855 }
856 
857 
858 void Foam::vtkPVFoam::PrintSelf(ostream& os, vtkIndent indent) const
859 {
860  os << indent << "Number of nodes: "
861  << (meshPtr_ ? meshPtr_->nPoints() : 0) << "\n";
862 
863  os << indent << "Number of cells: "
864  << (meshPtr_ ? meshPtr_->nCells() : 0) << "\n";
865 
866  os << indent << "Number of available time steps: "
867  << (dbPtr_.valid() ? dbPtr_().times().size() : 0) << "\n";
868 
869  os << indent << "mesh region: " << meshRegion_ << "\n";
870 }
871 
872 
873 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
string getEnv(const word &)
Return environment variable of given name.
Definition: POSIX.C:97
List< instant > instantList
List of instants.
Definition: instantList.H:42
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#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 void printMemory()
Simple memory used debugging information.
double * findTimes(int &nTimeSteps)
Allocate and return a list of selected times.
void updateInfo()
Update.
label nCells() const
void renderPatchNames(vtkRenderer *, const bool show)
Add/remove patch names to/from the view.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
void PrintSelf(ostream &, vtkIndent) const
Debug information.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
fileNameList findEtcFiles(const fileName &, bool mandatory=false, bool findFirst=false)
Search for files from user/group/shipped directories.
Definition: etcFiles.C:119
void Update(vtkMultiBlockDataSet *output, vtkMultiBlockDataSet *lagrangianOutput)
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
scalar y
void CleanUp()
Clean any storage.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:124
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
Functions to search &#39;etc&#39; directories for configuration files etc.
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
virtual fileName filePath(const bool checkGlobal, const IOobject &, const word &typeName) const =0
Search for an object. checkGlobal : also check undecomposed case.
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:205
const fileOperation & fileHandler()
Get current file handler.
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
static const char nl
Definition: Ostream.H:265
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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< string > stringList
A List of strings.
Definition: stringList.H:50
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
int setTime(int count, const double requestTimes[])
Set the runTime to the first plausible request time,.
fileName cwd()
Return current working directory path name.
Definition: POSIX.C:241
messageStream Info
vtkPVFoam(const char *const FileName, vtkPVFoamReader *reader)
Construct from components.
dimensioned< scalar > mag(const dimensioned< Type > &)
label nPoints() const
Field< vector > vectorField
Specialisation of Field<T> for vector.
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:253
~vtkPVFoam()
Destructor.
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
bool setEnv(const word &name, const std::string &value, const bool overwrite)
Set an environment variable.
Definition: POSIX.C:115
Namespace for OpenFOAM.
static label findClosestTimeIndex(const instantList &, const scalar, const word &constantName="constant")
Search instantList for the time index closest to the given time.
Definition: Time.C:750
IOporosityModelList pZones(mesh)