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