vtkPVblockMesh.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 "vtkPVblockMesh.H"
27 #include "vtkPVblockMeshReader.h"
28 
29 // OpenFOAM includes
30 #include "blockMesh.H"
31 #include "Time.H"
32 #include "patchZones.H"
33 #include "OStringStream.H"
34 
35 // VTK includes
36 #include "vtkDataArraySelection.h"
37 #include "vtkMultiBlockDataSet.h"
38 #include "vtkRenderer.h"
39 #include "vtkTextActor.h"
40 #include "vtkTextProperty.h"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  defineTypeNameAndDebug(vtkPVblockMesh, 0);
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::vtkPVblockMesh::resetCounters()
53 {
54  // Reset mesh part ids and sizes
55  arrayRangeBlocks_.reset();
56  arrayRangeEdges_.reset();
57  arrayRangeCorners_.reset();
58 }
59 
60 
61 void Foam::vtkPVblockMesh::updateInfoBlocks
62 (
63  vtkDataArraySelection* arraySelection
64 )
65 {
66  if (debug)
67  {
69  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "]" << endl;
70  }
71 
72  arrayRangeBlocks_.reset( arraySelection->GetNumberOfArrays() );
73 
74  const blockMesh& blkMesh = *meshPtr_;
75 
76  const int nBlocks = blkMesh.size();
77  for (int blockI = 0; blockI < nBlocks; ++blockI)
78  {
79  const blockDescriptor& blockDef = blkMesh[blockI];
80 
81  // Display either blockI as a number or with its name
82  // (looked up from blockMeshDict)
83  OStringStream os;
84  blockDescriptor::write(os, blockI, blkMesh.meshDict());
85  word partName(os.str());
86 
87  // append the (optional) zone name
88  if (!blockDef.zoneName().empty())
89  {
90  partName += " - " + blockDef.zoneName();
91  }
92 
93  // Add blockId and zoneName to GUI list
94  arraySelection->AddArray(partName.c_str());
95  }
96 
97  arrayRangeBlocks_ += nBlocks;
98 
99  if (debug)
100  {
101  // Just for debug info
102  getSelectedArrayEntries(arraySelection);
103  }
104 }
105 
106 
107 void Foam::vtkPVblockMesh::updateInfoEdges
108 (
109  vtkDataArraySelection* arraySelection
110 )
111 {
112  if (debug)
113  {
115  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "]" << endl;
116  }
117 
118  arrayRangeEdges_.reset( arraySelection->GetNumberOfArrays() );
119 
120  const blockMesh& blkMesh = *meshPtr_;
121  const blockEdgeList& edges = blkMesh.edges();
122 
123  const int nEdges = edges.size();
124  forAll(edges, edgeI)
125  {
126  OStringStream ostr;
127  blockVertex::write(ostr, edges[edgeI].start(), blkMesh.meshDict());
128  ostr<< ":";
129  blockVertex::write(ostr, edges[edgeI].end(), blkMesh.meshDict());
130  ostr << " - " << edges[edgeI].type();
131 
132  // Add "beg:end - type" to GUI list
133  arraySelection->AddArray(ostr.str().c_str());
134  }
135 
136  arrayRangeEdges_ += nEdges;
137 
138  if (debug)
139  {
140  // Just for debug info
141  getSelectedArrayEntries(arraySelection);
142  }
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
147 
149 (
150  const char* const FileName,
151  vtkPVblockMeshReader* reader
152 )
153 :
154  reader_(reader),
155  dbPtr_(nullptr),
156  meshPtr_(nullptr),
157  meshRegion_(polyMesh::defaultRegion),
158  meshDir_(polyMesh::meshSubDir),
159  arrayRangeBlocks_("block"),
160  arrayRangeEdges_("edges"),
161  arrayRangeCorners_("corners")
162 {
163  if (debug)
164  {
165  InfoInFunction<< " - " << FileName << endl;
166  }
167 
168  // Avoid argList and get rootPath/caseName directly from the file
169  fileName fullCasePath(fileName(FileName).path());
170 
171  if (!isDir(fullCasePath))
172  {
173  return;
174  }
175  if (fullCasePath == ".")
176  {
177  fullCasePath = cwd();
178  }
179 
180  // Set the case as an environment variable - some BCs might use this
181  if (fullCasePath.name().find("processor", 0) == 0)
182  {
183  const fileName globalCase = fullCasePath.path();
184 
185  setEnv("FOAM_CASE", globalCase, true);
186  setEnv("FOAM_CASENAME", globalCase.name(), true);
187  }
188  else
189  {
190  setEnv("FOAM_CASE", fullCasePath, true);
191  setEnv("FOAM_CASENAME", fullCasePath.name(), true);
192  }
193 
194  // Look for 'case{region}.OpenFOAM'
195  // could be stringent and insist the prefix match the directory name...
196  // Note: cannot use fileName::name() due to the embedded '{}'
197  string caseName(fileName(FileName).lessExt());
198  string::size_type beg = caseName.find_last_of("/{");
199  string::size_type end = caseName.find('}', beg);
200 
201  if
202  (
203  beg != string::npos && caseName[beg] == '{'
204  && end != string::npos && end == caseName.size()-1
205  )
206  {
207  meshRegion_ = caseName.substr(beg+1, end-beg-1);
208 
209  // Some safety
210  if (meshRegion_.empty())
211  {
212  meshRegion_ = polyMesh::defaultRegion;
213  }
214 
215  if (meshRegion_ != polyMesh::defaultRegion)
216  {
217  meshDir_ = meshRegion_/polyMesh::meshSubDir;
218  }
219  }
220 
221  if (debug)
222  {
223  Info<< " fullCasePath=" << fullCasePath << nl
224  << " FOAM_CASE=" << getEnv("FOAM_CASE") << nl
225  << " FOAM_CASENAME=" << getEnv("FOAM_CASENAME") << endl;
226  }
227 
228  // Create time object
229  dbPtr_.reset
230  (
231  new Time
232  (
233  Time::controlDictName,
234  fileName(fullCasePath.path()),
235  fileName(fullCasePath.name())
236  )
237  );
238 
239  dbPtr_().functionObjects().off();
240 
241  updateInfo();
242 }
243 
244 
245 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
246 
248 {
249  if (debug)
250  {
251  InfoInFunction << endl;
252  }
253 
254  // Hmm. pointNumberTextActors are not getting removed
255  //
256  forAll(pointNumberTextActorsPtrs_, pointi)
257  {
258  pointNumberTextActorsPtrs_[pointi]->Delete();
259  }
260  pointNumberTextActorsPtrs_.clear();
261 
262  delete meshPtr_;
263 }
264 
265 
266 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
267 
269 {
270  if (debug)
271  {
273  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "] " << endl;
274  }
275 
276  resetCounters();
277 
278  vtkDataArraySelection* blockSelection = reader_->GetBlockSelection();
279  vtkDataArraySelection* edgeSelection = reader_->GetCurvedEdgesSelection();
280 
281  // Enable 'internalMesh' on the first call
282  // or preserve the enabled selections
283  stringList enabledParts;
284  stringList enabledEdges;
285  bool firstTime = false;
286  if (!blockSelection->GetNumberOfArrays() && !meshPtr_)
287  {
288  firstTime = true;
289  }
290  else
291  {
292  enabledParts = getSelectedArrayEntries(blockSelection);
293  enabledEdges = getSelectedArrayEntries(edgeSelection);
294  }
295 
296  // Clear current mesh parts list
297  blockSelection->RemoveAllArrays();
298  edgeSelection->RemoveAllArrays();
299 
300  // Need a blockMesh
301  updateFoamMesh();
302 
303  // Update mesh parts list
304  updateInfoBlocks( blockSelection );
305 
306  // Update curved edges list
307  updateInfoEdges( edgeSelection );
308 
309  // Restore the enabled selections
310  if (!firstTime)
311  {
312  setSelectedArrayEntries(blockSelection, enabledParts);
313  setSelectedArrayEntries(edgeSelection, enabledEdges);
314  }
315 }
316 
317 
318 void Foam::vtkPVblockMesh::updateFoamMesh()
319 {
320  if (debug)
321  {
322  InfoInFunction << endl;
323  }
324 
325  // Check to see if the OpenFOAM mesh has been created
326  if (!meshPtr_)
327  {
328  if (debug)
329  {
331  << "Creating blockMesh at time=" << dbPtr_().timeName() << endl;
332  }
333 
334  // Set path for the blockMeshDict
335  const word dictName("blockMeshDict");
336  fileName dictPath(dbPtr_().system()/dictName);
337 
338  // Check if dictionary is present in the constant directory
339  if
340  (
341  exists
342  (
343  dbPtr_().path()/dbPtr_().constant()
345  )
346  )
347  {
348  dictPath = dbPtr_().constant()/polyMesh::meshSubDir/dictName;
349  }
350 
351  // Store dictionary since is used as database inside blockMesh class
352  // for names of vertices and blocks
353  IOdictionary* meshDictPtr = new IOdictionary
354  (
355  IOobject
356  (
357  dictPath,
358  dbPtr_(),
361  true
362  )
363  );
364  meshDictPtr->store();
365 
366  meshPtr_ = new blockMesh(*meshDictPtr, meshRegion_);
367  }
368 }
369 
370 
372 (
373  vtkMultiBlockDataSet* output
374 )
375 {
376  reader_->UpdateProgress(0.1);
377 
378  // Set up mesh parts selection(s)
379  updateBoolListStatus(blockStatus_, reader_->GetBlockSelection());
380 
381  // Set up curved edges selection(s)
382  updateBoolListStatus(edgeStatus_, reader_->GetCurvedEdgesSelection());
383 
384  reader_->UpdateProgress(0.2);
385 
386  // Update the OpenFOAM mesh
387  updateFoamMesh();
388  reader_->UpdateProgress(0.5);
389 
390  // Convert mesh element
391  int blockNo = 0;
392 
393  convertMeshCorners(output, blockNo);
394  convertMeshBlocks(output, blockNo);
395  convertMeshEdges(output, blockNo);
396 
397  reader_->UpdateProgress(0.8);
398 
399 }
400 
401 
403 {
404  reader_->UpdateProgress(1.0);
405 }
406 
407 
409 (
410  vtkRenderer* renderer,
411  const bool show
412 )
413 {
414  // Always remove old actors first
415 
416  forAll(pointNumberTextActorsPtrs_, pointi)
417  {
418  renderer->RemoveViewProp(pointNumberTextActorsPtrs_[pointi]);
419  pointNumberTextActorsPtrs_[pointi]->Delete();
420  }
421  pointNumberTextActorsPtrs_.clear();
422 
423  if (show && meshPtr_)
424  {
425  const blockMesh& blkMesh = *meshPtr_;
426  const pointField& cornerPts = blkMesh.vertices();
427  const scalar scaleFactor = blkMesh.scaleFactor();
428 
429  pointNumberTextActorsPtrs_.setSize(cornerPts.size());
430  forAll(cornerPts, pointi)
431  {
432  vtkTextActor* txt = vtkTextActor::New();
433 
434  // Display either pointi as a number or with its name
435  // (looked up from blockMeshDict)
436  {
437  OStringStream os;
438  blockVertex::write(os, pointi, blkMesh.meshDict());
439  txt->SetInput(os.str().c_str());
440  }
441 
442  // Set text properties
443  vtkTextProperty* tprop = txt->GetTextProperty();
444  tprop->SetFontFamilyToArial();
445  tprop->BoldOn();
446  tprop->ShadowOff();
447  tprop->SetLineSpacing(1.0);
448  tprop->SetFontSize(14);
449  tprop->SetColor(1.0, 0.0, 1.0);
450  tprop->SetJustificationToCentered();
451 
452  // Set text to use 3-D world co-ordinates
453  txt->GetPositionCoordinate()->SetCoordinateSystemToWorld();
454 
455  txt->GetPositionCoordinate()->SetValue
456  (
457  cornerPts[pointi].x()*scaleFactor,
458  cornerPts[pointi].y()*scaleFactor,
459  cornerPts[pointi].z()*scaleFactor
460  );
461 
462  // Add text to each renderer
463  renderer->AddViewProp(txt);
464 
465  // Maintain a list of text labels added so that they can be
466  // removed later
467  pointNumberTextActorsPtrs_[pointi] = txt;
468  }
469  }
470 }
471 
472 
473 
474 void Foam::vtkPVblockMesh::PrintSelf(ostream& os, vtkIndent indent) const
475 {
476 #if 0
477  os << indent << "Number of nodes: "
478  << (meshPtr_ ? meshPtr_->nPoints() : 0) << "\n";
479 
480  os << indent << "Number of cells: "
481  << (meshPtr_ ? meshPtr_->nCells() : 0) << "\n";
482 
483  os << indent << "Number of available time steps: "
484  << (dbPtr_.valid() ? dbPtr_().times().size() : 0) << endl;
485 #endif
486 }
487 
488 // ************************************************************************* //
string getEnv(const word &)
Return environment variable of given name.
Definition: POSIX.C:97
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
Definition: POSIX.C:520
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void PrintSelf(ostream &, vtkIndent) const
Debug information.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static void write(Ostream &, const label blocki, const dictionary &)
Write block index with dictionary lookup.
void renderPointNumbers(vtkRenderer *, const bool show)
Add/remove point numbers to/from the view.
scalar y
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
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
fileName dictPath
const word dictName("particleTrackDict")
graph_traits< Graph >::vertices_size_type size_type
Definition: SloanRenumber.C:73
void Update(vtkMultiBlockDataSet *output)
static const char nl
Definition: Ostream.H:260
defineTypeNameAndDebug(combustionModel, 0)
~vtkPVblockMesh()
Destructor.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
List< string > stringList
A List of strings.
Definition: stringList.H:50
void setSize(const label)
Reset size of List.
Definition: List.C:281
void CleanUp()
Clean any storage.
PtrList< blockEdge > blockEdgeList
A PtrList of blockEdges.
Definition: blockEdgeList.H:45
void updateInfo()
Update.
fileName cwd()
Return current working directory path name.
Definition: POSIX.C:241
vtkPVblockMesh(const char *const FileName, vtkPVblockMeshReader *reader)
Construct from components.
messageStream Info
static void write(Ostream &, const label, const dictionary &)
Write vertex index with optional name backsubstitution.
Definition: blockVertex.C:120
bool setEnv(const word &name, const std::string &value, const bool overwrite)
Set an environment variable.
Definition: POSIX.C:115
Namespace for OpenFOAM.
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1234
fileName path(UMean.rootPath()/UMean.caseName()/functionObjects::writeFile::outputPrefix/"graphs"/UMean.instance())
#define InfoInFunction
Report an information message using Foam::Info.