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