vtkPVblockMeshUtils.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 Description
25  Misc helper methods and utilities
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "vtkPVblockMesh.H"
30 #include "vtkPVblockMeshReader.h"
31 
32 // VTK includes
33 #include "vtkDataArraySelection.h"
34 #include "vtkDataSet.h"
35 #include "vtkMultiBlockDataSet.h"
36 #include "vtkInformation.h"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
43  // Extract up to the first non-word characters
44  inline word getFirstWord(const char* str)
45  {
46  if (str)
47  {
48  label n = 0;
49  while (str[n] && word::valid(str[n]))
50  {
51  ++n;
52  }
53  return word(str, n, true);
54  }
55  else
56  {
57  return word::null;
58  }
59 
60  }
62 
63 } // End namespace Foam
64 
65 
66 
67 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
68 
69 void Foam::vtkPVblockMesh::AddToBlock
70 (
71  vtkMultiBlockDataSet* output,
72  vtkDataSet* dataset,
73  const arrayRange& range,
74  const label datasetNo,
75  const std::string& datasetName
76 )
77 {
78  const int blockNo = range.block();
79 
80  vtkDataObject* blockDO = output->GetBlock(blockNo);
81  vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
82 
83  if (!block)
84  {
85  if (blockDO)
86  {
88  << "Block already has a vtkDataSet assigned to it"
89  << endl;
90  }
91 
92  block = vtkMultiBlockDataSet::New();
93  output->SetBlock(blockNo, block);
94  block->Delete();
95  }
96 
97  block->SetBlock(datasetNo, dataset);
98 
99  // name the block when assigning the first dataset
100  if (datasetNo == 0)
101  {
102  output->GetMetaData(blockNo)->Set
103  (
104  vtkCompositeDataSet::NAME(),
105  range.name()
106  );
107  }
108 
109  if (datasetName.size())
110  {
111  block->GetMetaData(datasetNo)->Set
112  (
113  vtkCompositeDataSet::NAME(),
114  datasetName.c_str()
115  );
116  }
117 }
118 
119 
120 vtkDataSet* Foam::vtkPVblockMesh::GetDataSetFromBlock
121 (
122  vtkMultiBlockDataSet* output,
123  const arrayRange& range,
124  const label datasetNo
125 )
126 {
127  const int blockNo = range.block();
128 
129  vtkDataObject* blockDO = output->GetBlock(blockNo);
130  vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
131 
132  if (block)
133  {
134  return vtkDataSet::SafeDownCast(block->GetBlock(datasetNo));
135  }
136 
137  return 0;
138 }
139 
140 
141 Foam::label Foam::vtkPVblockMesh::GetNumberOfDataSets
142 (
143  vtkMultiBlockDataSet* output,
144  const arrayRange& range
145 )
146 {
147  const int blockNo = range.block();
148 
149  vtkDataObject* blockDO = output->GetBlock(blockNo);
150  vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
151 
152  if (block)
153  {
154  return block->GetNumberOfBlocks();
155  }
156 
157  return 0;
158 }
159 
160 
161 Foam::wordHashSet Foam::vtkPVblockMesh::getSelected
162 (
163  vtkDataArraySelection* select
164 )
165 {
166  int nElem = select->GetNumberOfArrays();
167  wordHashSet selections(2*nElem);
168 
169  for (int elemI=0; elemI < nElem; ++elemI)
170  {
171  if (select->GetArraySetting(elemI))
172  {
173  selections.insert(getFirstWord(select->GetArrayName(elemI)));
174  }
175  }
176 
177  return selections;
178 }
179 
180 
181 Foam::wordHashSet Foam::vtkPVblockMesh::getSelected
182 (
183  vtkDataArraySelection* select,
184  const arrayRange& range
185 )
186 {
187  int nElem = select->GetNumberOfArrays();
188  wordHashSet selections(2*nElem);
189 
190  for (int elemI = range.start(); elemI < range.end(); ++elemI)
191  {
192  if (select->GetArraySetting(elemI))
193  {
194  selections.insert(getFirstWord(select->GetArrayName(elemI)));
195  }
196  }
197 
198  return selections;
199 }
200 
201 
202 Foam::stringList Foam::vtkPVblockMesh::getSelectedArrayEntries
203 (
204  vtkDataArraySelection* select,
205  const bool debug
206 )
207 {
208  stringList selections(select->GetNumberOfArrays());
209 
210  label nElem = 0;
211 
212  forAll(selections, elemI)
213  {
214  if (select->GetArraySetting(elemI))
215  {
216  selections[nElem++] = select->GetArrayName(elemI);
217  }
218  }
219 
220  selections.setSize(nElem);
221 
222  // Debugging ...
223  const label nAvailableElem = select->GetNumberOfArrays();
225  << "available(";
226  for (int elemI = 0; elemI < nAvailableElem; ++elemI)
227  {
228  DebugInfo<< " \"" << select->GetArrayName(elemI) << "\"";
229  }
230  DebugInfo<< " )\nselected(";
231  forAll(selections, elemI)
232  {
233  DebugInfo<< " " << selections[elemI];
234  }
235  DebugInfo<< " )" << endl;
236 
237  return selections;
238 }
239 
240 
241 Foam::stringList Foam::vtkPVblockMesh::getSelectedArrayEntries
242 (
243  vtkDataArraySelection* select,
244  const arrayRange& range,
245  const bool debug
246 )
247 {
248  stringList selections(range.size());
249 
250  label nElem = 0;
251 
252  for (int elemI = range.start(); elemI < range.end(); ++elemI)
253  {
254  if (select->GetArraySetting(elemI))
255  {
256  selections[nElem++] = select->GetArrayName(elemI);
257  }
258  }
259 
260  selections.setSize(nElem);
261 
262  // Debugging ...
264  << "available(";
265  for (int elemI = range.start(); elemI < range.end(); ++elemI)
266  {
267  DebugInfo<< " \"" << select->GetArrayName(elemI) << "\"";
268  }
269  DebugInfo<< " )\nselected(";
270  forAll(selections, elemI)
271  {
272  DebugInfo<< " " << selections[elemI];
273  }
274  DebugInfo<< " )" << endl;
275 
276  return selections;
277 }
278 
279 
280 void Foam::vtkPVblockMesh::setSelectedArrayEntries
281 (
282  vtkDataArraySelection* select,
283  const stringList& selections
284 )
285 {
286  const int nElem = select->GetNumberOfArrays();
287  select->DisableAllArrays();
288 
289  // Loop through entries, setting values from selectedEntries
290  for (int elemI=0; elemI < nElem; ++elemI)
291  {
292  string arrayName(select->GetArrayName(elemI));
293 
294  forAll(selections, elemI)
295  {
296  if (selections[elemI] == arrayName)
297  {
298  select->EnableArray(arrayName.c_str());
299  break;
300  }
301  }
302  }
303 }
304 
305 
306 void Foam::vtkPVblockMesh::updateBoolListStatus
307 (
308  boolList& status,
309  vtkDataArraySelection* selection
310 )
311 {
312  const label nElem = selection->GetNumberOfArrays();
313 
314  if (status.size() != nElem)
315  {
316  status.setSize(nElem);
317  status = false;
318  }
319 
320  forAll(status, elemI)
321  {
322  status[elemI] = selection->GetArraySetting(elemI);
323  }
324 }
325 
326 
327 // ************************************************************************* //
scalar range
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
A HashTable with keys but without contents.
Definition: HashSet.H:62
static bool valid(char)
Is this character valid for a word.
Definition: wordI.H:115
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
#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.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
List< string > stringList
A List of strings.
Definition: stringList.H:50