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