All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 
220 Foam::vtkPVFoam::vtkPVFoam
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  // Create time object
326  dbPtr_.reset
327  (
328  new Time
329  (
330  Time::controlDictName,
331  fileName(fullCasePath.path()),
332  fileName(fullCasePath.name())
333  )
334  );
335 
336  dbPtr_().functionObjects().off();
337 
338  fileNameList configDictFiles = findEtcFiles("paraFoam", false);
339  forAllReverse(configDictFiles, cdfi)
340  {
341  configDict_.merge(dictionary(IFstream(configDictFiles[cdfi])()));
342  }
343 
344  updateInfo();
345 }
346 
347 
348 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
349 
351 {
352  if (debug)
353  {
354  Info<< "<end> Foam::vtkPVFoam::~vtkPVFoam" << endl;
355  }
356 
357  delete meshPtr_;
358 }
359 
360 
361 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
362 
364 {
365  if (debug)
366  {
367  Info<< "<beg> Foam::vtkPVFoam::updateInfo"
368  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "] timeIndex="
369  << timeIndex_ << endl;
370  }
371 
372  resetCounters();
373 
374  vtkDataArraySelection* partSelection = reader_->GetPartSelection();
375 
376  // there are two ways to ensure we have the correct list of parts:
377  // 1. remove everything and then set particular entries 'on'
378  // 2. build a 'char **' list and call SetArraysWithDefault()
379  //
380  // Nr. 2 has the potential advantage of not touching the modification
381  // time of the vtkDataArraySelection, but the qt/paraview proxy
382  // layer doesn't care about that anyhow.
383 
384  // enable 'internalMesh' on the first call
385  // or preserve the enabled selections
386  stringList enabledEntries;
387  bool first = !partSelection->GetNumberOfArrays() && !meshPtr_;
388  if (first)
389  {
390  enabledEntries.setSize(1);
391  enabledEntries[0] = "internalMesh";
392  }
393  else
394  {
395  enabledEntries = getSelectedArrayEntries(partSelection);
396  }
397 
398  // Clear current mesh parts list
399  partSelection->RemoveAllArrays();
400 
401  // Update mesh parts list - add Lagrangian at the bottom
402  updateInfoInternalMesh(partSelection);
403  updateInfoPatches(partSelection, enabledEntries, first);
404  updateInfoSets(partSelection);
405  updateInfoZones(partSelection);
406  updateInfoLagrangian(partSelection);
407 
408  // restore the enabled selections
409  setSelectedArrayEntries(partSelection, enabledEntries);
410 
411  if (meshChanged_)
412  {
413  fieldsChanged_ = true;
414  }
415 
416  // Update volume, point and lagrangian fields
417  updateInfoFields<fvPatchField, volMesh>
418  (
419  reader_->GetVolFieldSelection()
420  );
421  updateInfoFields<pointPatchField, pointMesh>
422  (
423  reader_->GetPointFieldSelection()
424  );
425  updateInfoLagrangianFields();
426 
427  if (debug)
428  {
429  // just for debug info
430  getSelectedArrayEntries(partSelection);
431  Info<< "<end> Foam::vtkPVFoam::updateInfo" << endl;
432  }
433 
434 }
435 
436 
437 void Foam::vtkPVFoam::updateFoamMesh()
438 {
439  if (debug)
440  {
441  Info<< "<beg> Foam::vtkPVFoam::updateFoamMesh" << endl;
442  printMemory();
443  }
444 
445  if (!reader_->GetCacheMesh())
446  {
447  delete meshPtr_;
448  meshPtr_ = nullptr;
449  }
450 
451  // Check to see if the OpenFOAM mesh has been created
452  if (!meshPtr_)
453  {
454  if (debug)
455  {
456  Info<< "Creating OpenFOAM mesh for region " << meshRegion_
457  << " at time=" << dbPtr_().timeName()
458  << endl;
459 
460  }
461 
462  meshPtr_ = new fvMesh
463  (
464  IOobject
465  (
466  meshRegion_,
467  dbPtr_().timeName(),
468  dbPtr_(),
470  )
471  );
472 
473  meshChanged_ = true;
474  }
475  else
476  {
477  if (debug)
478  {
479  Info<< "Using existing OpenFOAM mesh" << endl;
480  }
481  }
482 
483  if (debug)
484  {
485  Info<< "<end> Foam::vtkPVFoam::updateFoamMesh" << endl;
486  printMemory();
487  }
488 }
489 
490 
492 (
493  vtkMultiBlockDataSet* output,
494  vtkMultiBlockDataSet* lagrangianOutput
495 )
496 {
497  if (debug)
498  {
499  cout<< "<beg> Foam::vtkPVFoam::Update - output with "
500  << output->GetNumberOfBlocks() << " and "
501  << lagrangianOutput->GetNumberOfBlocks() << " blocks\n";
502  output->Print(cout);
503  lagrangianOutput->Print(cout);
504  printMemory();
505  }
506  reader_->UpdateProgress(0.1);
507 
508  // Set up mesh parts selection(s)
509  updateMeshPartsStatus();
510 
511  reader_->UpdateProgress(0.15);
512 
513  // Update the OpenFOAM mesh
514  updateFoamMesh();
515  reader_->UpdateProgress(0.4);
516 
517  // Convert meshes - start port0 at block=0
518  int blockNo = 0;
519 
520  convertMeshVolume(output, blockNo);
521  convertMeshPatches(output, blockNo);
522  reader_->UpdateProgress(0.6);
523 
524  if (reader_->GetIncludeZones())
525  {
526  convertMeshCellZones(output, blockNo);
527  convertMeshFaceZones(output, blockNo);
528  convertMeshPointZones(output, blockNo);
529  reader_->UpdateProgress(0.65);
530  }
531 
532  if (reader_->GetIncludeSets())
533  {
534  convertMeshCellSets(output, blockNo);
535  convertMeshFaceSets(output, blockNo);
536  convertMeshPointSets(output, blockNo);
537  reader_->UpdateProgress(0.7);
538  }
539 
540  convertMeshLagrangian(lagrangianOutput, blockNo);
541 
542  reader_->UpdateProgress(0.8);
543 
544  // Update fields
545  convertVolFields(output);
546  convertPointFields(output);
547  convertLagrangianFields(lagrangianOutput);
548  if (debug)
549  {
550  Info<< "done reader part" << endl;
551  }
552  reader_->UpdateProgress(0.95);
553 
554  meshChanged_ = fieldsChanged_ = false;
555 }
556 
557 
559 {
560  // reclaim some memory
561  reduceMemory();
562  reader_->UpdateProgress(1.0);
563 }
564 
565 
566 double* Foam::vtkPVFoam::findTimes(int& nTimeSteps)
567 {
568  int nTimes = 0;
569  double* tsteps = nullptr;
570 
571  if (dbPtr_.valid())
572  {
573  Time& runTime = dbPtr_();
574  // Get times list. Flush first to force refresh.
575  fileHandler().flush();
576  instantList timeLst = runTime.times();
577 
578  // find the first time for which this mesh appears to exist
579  label timeI = 0;
580  for (; timeI < timeLst.size(); ++timeI)
581  {
582  const word& timeName = timeLst[timeI].name();
583 
584  if
585  (
586  IOobject
587  (
588  "points",
589  timeName,
590  meshDir_,
591  runTime
592  ).typeHeaderOk<pointIOField>(true)
593  )
594  {
595  break;
596  }
597  }
598 
599  nTimes = timeLst.size() - timeI;
600 
601  // skip "constant" time whenever possible
602  if (timeI == 0 && nTimes > 1)
603  {
604  if (timeLst[timeI].name() == runTime.constant())
605  {
606  ++timeI;
607  --nTimes;
608  }
609  }
610 
611 
612  // skip "0/" time if requested and possible
613  if (nTimes > 1 && reader_->GetSkipZeroTime())
614  {
615  if (mag(timeLst[timeI].value()) < small)
616  {
617  ++timeI;
618  --nTimes;
619  }
620  }
621 
622  if (nTimes)
623  {
624  tsteps = new double[nTimes];
625  for (label stepI = 0; stepI < nTimes; ++stepI, ++timeI)
626  {
627  tsteps[stepI] = timeLst[timeI].value();
628  }
629  }
630  }
631  else
632  {
633  if (debug)
634  {
635  cout<< "no valid dbPtr:\n";
636  }
637  }
638 
639  // vector length returned via the parameter
640  nTimeSteps = nTimes;
641 
642  return tsteps;
643 }
644 
645 
647 (
648  vtkRenderer* renderer,
649  const bool show
650 )
651 {
652  if (!meshPtr_)
653  {
654  return;
655  }
656 
657  // always remove old actors first
658 
659  forAll(patchTextActorsPtrs_, patchi)
660  {
661  renderer->RemoveViewProp(patchTextActorsPtrs_[patchi]);
662  patchTextActorsPtrs_[patchi]->Delete();
663  }
664  patchTextActorsPtrs_.clear();
665 
666  if (show)
667  {
668  // get the display patches, strip off any suffix
669  wordHashSet selectedPatches = getSelected
670  (
671  reader_->GetPartSelection(),
672  arrayRangePatches_
673  );
674 
675  if (selectedPatches.empty())
676  {
677  return;
678  }
679 
680  const polyBoundaryMesh& pbMesh = meshPtr_->boundaryMesh();
681 
682  // Find the total number of zones
683  // Each zone will take the patch name
684  // Number of zones per patch ... zero zones should be skipped
685  labelList nZones(pbMesh.size(), 0);
686 
687  // Per global zone number the average face centre position
688  List<DynamicList<point>> zoneCentre(pbMesh.size());
689 
690 
691  // Loop through all patches to determine zones, and centre of each zone
692  forAll(pbMesh, patchi)
693  {
694  const polyPatch& pp = pbMesh[patchi];
695 
696  // Only include the patch if it is selected
697  if (!selectedPatches.found(pp.name()))
698  {
699  continue;
700  }
701 
702  const labelListList& edgeFaces = pp.edgeFaces();
703  const vectorField& n = pp.faceNormals();
704 
705  boolList featEdge(pp.nEdges(), false);
706 
707  forAll(edgeFaces, edgeI)
708  {
709  const labelList& eFaces = edgeFaces[edgeI];
710 
711  if (eFaces.size() == 1)
712  {
713  // Note: could also do ones with > 2 faces but this gives
714  // too many zones for baffles
715  featEdge[edgeI] = true;
716  }
717  else if (mag(n[eFaces[0]] & n[eFaces[1]]) < 0.5)
718  {
719  featEdge[edgeI] = true;
720  }
721  }
722 
723  // Do topological analysis of patch, find disconnected regions
724  patchZones pZones(pp, featEdge);
725 
726  nZones[patchi] = pZones.nZones();
727 
728  labelList zoneNFaces(pZones.nZones(), 0);
729 
730  // Create storage for additional zone centres
731  forAll(zoneNFaces, zoneI)
732  {
733  zoneCentre[patchi].append(Zero);
734  }
735 
736  // Do averaging per individual zone
737  forAll(pp, facei)
738  {
739  label zoneI = pZones[facei];
740  zoneCentre[patchi][zoneI] += pp[facei].centre(pp.points());
741  zoneNFaces[zoneI]++;
742  }
743 
744  forAll(zoneCentre[patchi], zoneI)
745  {
746  zoneCentre[patchi][zoneI] /= zoneNFaces[zoneI];
747  }
748  }
749 
750  // Count number of zones we're actually going to display.
751  // This is truncated to a max per patch
752 
753  const label MAXPATCHZONES = 20;
754 
755  label displayZoneI = 0;
756 
757  forAll(pbMesh, patchi)
758  {
759  displayZoneI += min(MAXPATCHZONES, nZones[patchi]);
760  }
761 
762  if (debug)
763  {
764  Info<< "displayed zone centres = " << displayZoneI << nl
765  << "zones per patch = " << nZones << endl;
766  }
767 
768  // Set the size of the patch labels to max number of zones
769  patchTextActorsPtrs_.setSize(displayZoneI);
770 
771  if (debug)
772  {
773  Info<< "constructing patch labels" << endl;
774  }
775 
776  // Actor index
777  displayZoneI = 0;
778 
779  forAll(pbMesh, patchi)
780  {
781  const polyPatch& pp = pbMesh[patchi];
782 
783  label globalZoneI = 0;
784 
785  // Only selected patches will have a non-zero number of zones
786  label nDisplayZones = min(MAXPATCHZONES, nZones[patchi]);
787  label increment = 1;
788  if (nZones[patchi] >= MAXPATCHZONES)
789  {
790  increment = nZones[patchi]/MAXPATCHZONES;
791  }
792 
793  for (label i = 0; i < nDisplayZones; i++)
794  {
795  if (debug)
796  {
797  Info<< "patch name = " << pp.name() << nl
798  << "anchor = " << zoneCentre[patchi][globalZoneI] << nl
799  << "globalZoneI = " << globalZoneI << endl;
800  }
801 
802  vtkTextActor* txt = vtkTextActor::New();
803 
804  txt->SetInput(pp.name().c_str());
805 
806  // Set text properties
807  vtkTextProperty* tprop = txt->GetTextProperty();
808  tprop->SetFontFamilyToArial();
809  tprop->BoldOff();
810  tprop->ShadowOff();
811  tprop->SetLineSpacing(1.0);
812  tprop->SetFontSize(12);
813  tprop->SetColor(1.0, 0.0, 0.0);
814  tprop->SetJustificationToCentered();
815 
816  // Set text to use 3-D world co-ordinates
817  txt->GetPositionCoordinate()->SetCoordinateSystemToWorld();
818 
819  txt->GetPositionCoordinate()->SetValue
820  (
821  zoneCentre[patchi][globalZoneI].x(),
822  zoneCentre[patchi][globalZoneI].y(),
823  zoneCentre[patchi][globalZoneI].z()
824  );
825 
826  // Add text to each renderer
827  renderer->AddViewProp(txt);
828 
829  // Maintain a list of text labels added so that they can be
830  // removed later
831  patchTextActorsPtrs_[displayZoneI] = txt;
832 
833  globalZoneI += increment;
834  displayZoneI++;
835  }
836  }
837 
838  // Resize the patch names list to the actual number of patch names added
839  patchTextActorsPtrs_.setSize(displayZoneI);
840  }
841 }
842 
843 
844 void Foam::vtkPVFoam::PrintSelf(ostream& os, vtkIndent indent) const
845 {
846  os << indent << "Number of nodes: "
847  << (meshPtr_ ? meshPtr_->nPoints() : 0) << "\n";
848 
849  os << indent << "Number of cells: "
850  << (meshPtr_ ? meshPtr_->nCells() : 0) << "\n";
851 
852  os << indent << "Number of available time steps: "
853  << (dbPtr_.valid() ? dbPtr_().times().size() : 0) << "\n";
854 
855  os << indent << "mesh region: " << meshRegion_ << "\n";
856 }
857 
858 
859 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
string getEnv(const word &)
Return environment variable of given name.
Definition: POSIX.C:96
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: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 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:440
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:528
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:207
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:240
messageStream Info
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:249
~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:114
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)