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