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-2021 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 "OSspecific.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  {
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 }
106 
107 
108 void Foam::vtkPVblockMesh::updateInfoEdges
109 (
110  vtkDataArraySelection* arraySelection
111 )
112 {
113  if (debug)
114  {
116  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "]" << endl;
117  }
118 
119  arrayRangeEdges_.reset( arraySelection->GetNumberOfArrays() );
120 
121  const blockMesh& blkMesh = *meshPtr_;
122  const blockEdgeList& edges = blkMesh.edges();
123 
124  const int nEdges = edges.size();
125  forAll(edges, edgeI)
126  {
127  OStringStream ostr;
128  blockVertex::write(ostr, edges[edgeI].start(), blkMesh.meshDict());
129  ostr<< ":";
130  blockVertex::write(ostr, edges[edgeI].end(), blkMesh.meshDict());
131  ostr << " - " << edges[edgeI].type();
132 
133  // Add "beg:end - type" to GUI list
134  arraySelection->AddArray(ostr.str().c_str());
135  }
136 
137  arrayRangeEdges_ += nEdges;
138 
139  if (debug)
140  {
141  // Just for debug info
142  getSelectedArrayEntries(arraySelection);
143  }
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
148 
150 (
151  const char* const FileName,
152  vtkPVblockMeshReader* reader
153 )
154 :
155  reader_(reader),
156  dbPtr_(nullptr),
157  meshPtr_(nullptr),
158  meshRegion_(polyMesh::defaultRegion),
159  meshDir_(polyMesh::meshSubDir),
160  arrayRangeBlocks_("block"),
161  arrayRangeEdges_("edges"),
162  arrayRangeCorners_("corners")
163 {
164  if (debug)
165  {
166  InfoInFunction<< " - " << FileName << endl;
167  }
168 
169  // Avoid argList and get rootPath/caseName directly from the file
170  fileName fullCasePath(fileName(FileName).path());
171 
172  if (!isDir(fullCasePath))
173  {
174  return;
175  }
176  if (fullCasePath == ".")
177  {
178  fullCasePath = cwd();
179  }
180 
181  // Set the case as an environment variable - some BCs might use this
182  if (fullCasePath.name().find("processor", 0) == 0)
183  {
184  const fileName globalCase = fullCasePath.path();
185 
186  setEnv("FOAM_CASE", globalCase, true);
187  setEnv("FOAM_CASENAME", globalCase.name(), true);
188  }
189  else
190  {
191  setEnv("FOAM_CASE", fullCasePath, true);
192  setEnv("FOAM_CASENAME", fullCasePath.name(), true);
193  }
194 
195  // Look for 'case{region}.OpenFOAM'
196  // could be stringent and insist the prefix match the directory name...
197  // Note: cannot use fileName::name() due to the embedded '{}'
198  string caseName(fileName(FileName).lessExt());
199  string::size_type beg = caseName.find_last_of("/{");
200  string::size_type end = caseName.find('}', beg);
201 
202  if
203  (
204  beg != string::npos && caseName[beg] == '{'
205  && end != string::npos && end == caseName.size()-1
206  )
207  {
208  meshRegion_ = caseName.substr(beg+1, end-beg-1);
209 
210  // Some safety
211  if (meshRegion_.empty())
212  {
213  meshRegion_ = polyMesh::defaultRegion;
214  }
215 
216  if (meshRegion_ != polyMesh::defaultRegion)
217  {
218  meshDir_ = meshRegion_/polyMesh::meshSubDir;
219  }
220  }
221 
222  if (debug)
223  {
224  Info<< " fullCasePath=" << fullCasePath << nl
225  << " FOAM_CASE=" << getEnv("FOAM_CASE") << nl
226  << " FOAM_CASENAME=" << getEnv("FOAM_CASENAME") << endl;
227  }
228 
229  // Create time object
230  dbPtr_.reset
231  (
232  new Time
233  (
234  Time::controlDictName,
235  fileName(fullCasePath.path()),
236  fileName(fullCasePath.name())
237  )
238  );
239 
240  dbPtr_().functionObjects().off();
241 
242  updateInfo();
243 }
244 
245 
246 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
247 
249 {
250  if (debug)
251  {
252  InfoInFunction << endl;
253  }
254 
255  // Hmm. pointNumberTextActors are not getting removed
256  //
257  forAll(pointNumberTextActorsPtrs_, pointi)
258  {
259  pointNumberTextActorsPtrs_[pointi]->Delete();
260  }
261  pointNumberTextActorsPtrs_.clear();
262 
263  delete meshPtr_;
264 }
265 
266 
267 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
268 
270 {
271  if (debug)
272  {
274  << " [meshPtr=" << (meshPtr_ ? "set" : "nullptr") << "] " << endl;
275  }
276 
277  resetCounters();
278 
279  vtkDataArraySelection* blockSelection = reader_->GetBlockSelection();
280  vtkDataArraySelection* edgeSelection = reader_->GetCurvedEdgesSelection();
281 
282  // Enable 'internalMesh' on the first call
283  // or preserve the enabled selections
284  stringList enabledParts;
285  stringList enabledEdges;
286  bool firstTime = false;
287  if (!blockSelection->GetNumberOfArrays() && !meshPtr_)
288  {
289  firstTime = true;
290  }
291  else
292  {
293  enabledParts = getSelectedArrayEntries(blockSelection);
294  enabledEdges = getSelectedArrayEntries(edgeSelection);
295  }
296 
297  // Clear current mesh parts list
298  blockSelection->RemoveAllArrays();
299  edgeSelection->RemoveAllArrays();
300 
301  // Need a blockMesh
302  updateFoamMesh();
303 
304  // Update mesh parts list
305  updateInfoBlocks( blockSelection );
306 
307  // Update curved edges list
308  updateInfoEdges( edgeSelection );
309 
310  // Restore the enabled selections
311  if (!firstTime)
312  {
313  setSelectedArrayEntries(blockSelection, enabledParts);
314  setSelectedArrayEntries(edgeSelection, enabledEdges);
315  }
316 }
317 
318 
319 void Foam::vtkPVblockMesh::updateFoamMesh()
320 {
321  if (debug)
322  {
323  InfoInFunction << endl;
324  }
325 
326  // Check to see if the OpenFOAM mesh has been created
327  if (!meshPtr_)
328  {
329  if (debug)
330  {
332  << "Creating blockMesh at time=" << dbPtr_().timeName() << endl;
333  }
334 
335  // Set path for the blockMeshDict
336  const word dictName("blockMeshDict");
337  fileName dictPath(dbPtr_().system()/dictName);
338 
339  // Check if dictionary is present in the constant directory
340  if
341  (
342  exists
343  (
344  dbPtr_().path()/dbPtr_().constant()
346  )
347  )
348  {
349  dictPath = dbPtr_().constant()/polyMesh::meshSubDir/dictName;
350  }
351 
352  // Store dictionary since is used as database inside blockMesh class
353  // for names of vertices and blocks
354  IOdictionary* meshDictPtr = new IOdictionary
355  (
356  IOobject
357  (
358  dictPath,
359  dbPtr_(),
362  true
363  )
364  );
365  meshDictPtr->store();
366 
367  meshPtr_ = new blockMesh(*meshDictPtr, meshRegion_);
368  }
369 }
370 
371 
373 (
374  vtkMultiBlockDataSet* output
375 )
376 {
377  reader_->UpdateProgress(0.1);
378 
379  // Set up mesh parts selection(s)
380  updateBoolListStatus(blockStatus_, reader_->GetBlockSelection());
381 
382  // Set up curved edges selection(s)
383  updateBoolListStatus(edgeStatus_, reader_->GetCurvedEdgesSelection());
384 
385  reader_->UpdateProgress(0.2);
386 
387  // Update the OpenFOAM mesh
388  updateFoamMesh();
389  reader_->UpdateProgress(0.5);
390 
391  // Convert mesh element
392  int blockNo = 0;
393 
394  convertMeshCorners(output, blockNo);
395  convertMeshBlocks(output, blockNo);
396  convertMeshEdges(output, blockNo);
397 
398  reader_->UpdateProgress(0.8);
399 
400 }
401 
402 
404 {
405  reader_->UpdateProgress(1.0);
406 }
407 
408 
410 (
411  vtkRenderer* renderer,
412  const bool show
413 )
414 {
415  // Always remove old actors first
416 
417  forAll(pointNumberTextActorsPtrs_, pointi)
418  {
419  renderer->RemoveViewProp(pointNumberTextActorsPtrs_[pointi]);
420  pointNumberTextActorsPtrs_[pointi]->Delete();
421  }
422  pointNumberTextActorsPtrs_.clear();
423 
424  if (show && meshPtr_)
425  {
426  const blockMesh& blkMesh = *meshPtr_;
427  const pointField& cornerPts = blkMesh.vertices();
428  const scalar scaleFactor = blkMesh.scaleFactor();
429 
430  pointNumberTextActorsPtrs_.setSize(cornerPts.size());
431  forAll(cornerPts, pointi)
432  {
433  vtkTextActor* txt = vtkTextActor::New();
434 
435  // Display either pointi as a number or with its name
436  // (looked up from blockMeshDict)
437  {
438  OStringStream os;
439  blockVertex::write(os, pointi, blkMesh.meshDict());
440  txt->SetInput(os.str().c_str());
441  }
442 
443  // Set text properties
444  vtkTextProperty* tprop = txt->GetTextProperty();
445  tprop->SetFontFamilyToArial();
446  tprop->BoldOn();
447  tprop->ShadowOff();
448  tprop->SetLineSpacing(1.0);
449  tprop->SetFontSize(14);
450  tprop->SetColor(1.0, 0.0, 1.0);
451  tprop->SetJustificationToCentered();
452 
453  // Set text to use 3-D world co-ordinates
454  txt->GetPositionCoordinate()->SetCoordinateSystemToWorld();
455 
456  txt->GetPositionCoordinate()->SetValue
457  (
458  cornerPts[pointi].x()*scaleFactor,
459  cornerPts[pointi].y()*scaleFactor,
460  cornerPts[pointi].z()*scaleFactor
461  );
462 
463  // Add text to each renderer
464  renderer->AddViewProp(txt);
465 
466  // Maintain a list of text labels added so that they can be
467  // removed later
468  pointNumberTextActorsPtrs_[pointi] = txt;
469  }
470  }
471 }
472 
473 
474 
475 void Foam::vtkPVblockMesh::PrintSelf(ostream& os, vtkIndent indent) const
476 {
477 #if 0
478  os << indent << "Number of nodes: "
479  << (meshPtr_ ? meshPtr_->nPoints() : 0) << "\n";
480 
481  os << indent << "Number of cells: "
482  << (meshPtr_ ? meshPtr_->nCells() : 0) << "\n";
483 
484  os << indent << "Number of available time steps: "
485  << (dbPtr_.valid() ? dbPtr_().times().size() : 0) << endl;
486 #endif
487 }
488 
489 // ************************************************************************* //
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
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
void PrintSelf(ostream &, vtkIndent) const
Debug information.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
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.
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
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:125
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
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
const word dictName("noiseDict")
Namespace for OpenFOAM.
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1230
#define InfoInFunction
Report an information message using Foam::Info.