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