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