vtkPVFoamUtils.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 "vtkPVFoam.H"
30 #include "vtkPVFoamReader.h"
31 
32 // OpenFOAM includes
33 #include "fvMesh.H"
34 #include "Time.H"
35 #include "IFstream.H"
36 
37 // VTK includes
38 #include "vtkDataArraySelection.h"
39 #include "vtkDataSet.h"
40 #include "vtkMultiBlockDataSet.h"
41 #include "vtkInformation.h"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47  // Extract up to the first non-word characters
48  inline word getFirstWord(const char* str)
49  {
50  if (str)
51  {
52  label n = 0;
53  while (str[n] && word::valid(str[n]))
54  {
55  ++n;
56  }
57  return word(str, n, true);
58  }
59  else
60  {
61  return word::null;
62  }
63 
64  }
65 } // End namespace Foam
66 
67 
68 
69 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
70 
71 void Foam::vtkPVFoam::AddToBlock
72 (
73  vtkMultiBlockDataSet* output,
74  vtkDataSet* dataset,
75  const arrayRange& range,
76  const label datasetNo,
77  const std::string& datasetName
78 )
79 {
80  const int blockNo = range.block();
81 
82  vtkDataObject* blockDO = output->GetBlock(blockNo);
83  vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
84 
85  if (!block)
86  {
87  if (blockDO)
88  {
90  << "Block already has a vtkDataSet assigned to it"
91  << endl;
92  }
93 
94  block = vtkMultiBlockDataSet::New();
95  output->SetBlock(blockNo, block);
96  block->Delete();
97  }
98 
99  block->SetBlock(datasetNo, dataset);
100 
101  // name the block when assigning the first dataset
102  if (datasetNo == 0)
103  {
104  output->GetMetaData(blockNo)->Set
105  (
106  vtkCompositeDataSet::NAME(),
107  range.name()
108  );
109  }
110 
111  if (datasetName.size())
112  {
113  block->GetMetaData(datasetNo)->Set
114  (
115  vtkCompositeDataSet::NAME(),
116  datasetName.c_str()
117  );
118  }
119 }
120 
121 
122 vtkDataSet* Foam::vtkPVFoam::GetDataSetFromBlock
123 (
124  vtkMultiBlockDataSet* output,
125  const arrayRange& range,
126  const label datasetNo
127 )
128 {
129  const int blockNo = range.block();
130 
131  vtkDataObject* blockDO = output->GetBlock(blockNo);
132  vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
133 
134  if (block)
135  {
136  return vtkDataSet::SafeDownCast(block->GetBlock(datasetNo));
137  }
138 
139  return 0;
140 }
141 
142 
143 Foam::label Foam::vtkPVFoam::GetNumberOfDataSets
144 (
145  vtkMultiBlockDataSet* output,
146  const arrayRange& range
147 )
148 {
149  const int blockNo = range.block();
150 
151  vtkDataObject* blockDO = output->GetBlock(blockNo);
152  vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
153  if (block)
154  {
155  return block->GetNumberOfBlocks();
156  }
157 
158  return 0;
159 }
160 
161 
162 Foam::word Foam::vtkPVFoam::getPartName(const int partId)
163 {
164  return getFirstWord(reader_->GetPartArrayName(partId));
165 }
166 
167 
168 Foam::wordHashSet Foam::vtkPVFoam::getSelected
169 (
170  vtkDataArraySelection* select
171 )
172 {
173  int nElem = select->GetNumberOfArrays();
174  wordHashSet selections(2*nElem);
175 
176  for (int elemI=0; elemI < nElem; ++elemI)
177  {
178  if (select->GetArraySetting(elemI))
179  {
180  selections.insert(getFirstWord(select->GetArrayName(elemI)));
181  }
182  }
183 
184  return selections;
185 }
186 
187 
188 Foam::wordHashSet Foam::vtkPVFoam::getSelected
189 (
190  vtkDataArraySelection* select,
191  const arrayRange& range
192 )
193 {
194  int nElem = select->GetNumberOfArrays();
195  wordHashSet selections(2*nElem);
196 
197  for (int elemI = range.start(); elemI < range.end(); ++elemI)
198  {
199  if (select->GetArraySetting(elemI))
200  {
201  selections.insert(getFirstWord(select->GetArrayName(elemI)));
202  }
203  }
204 
205  return selections;
206 }
207 
208 
209 Foam::stringList Foam::vtkPVFoam::getSelectedArrayEntries
210 (
211  vtkDataArraySelection* select,
212  const bool debug
213 )
214 {
215  stringList selections(select->GetNumberOfArrays());
216 
217  label nElem = 0;
218 
219  forAll(selections, elemI)
220  {
221  if (select->GetArraySetting(elemI))
222  {
223  selections[nElem++] = select->GetArrayName(elemI);
224  }
225  }
226 
227  selections.setSize(nElem);
228 
229  // Debugging ...
230  const label nAvailableElem = select->GetNumberOfArrays();
231  DebugInFunction<< "available(";
232  for (int elemI = 0; elemI < nAvailableElem; ++elemI)
233  {
234  DebugInfo<< " \"" << select->GetArrayName(elemI) << "\"";
235  }
236  DebugInfo<< " )\n selected(";
237  forAll(selections, elemI)
238  {
239  DebugInfo<< " " << selections[elemI];
240  }
241  DebugInfo<< " )" << endl;
242 
243  return selections;
244 }
245 
246 
247 Foam::stringList Foam::vtkPVFoam::getSelectedArrayEntries
248 (
249  vtkDataArraySelection* select,
250  const arrayRange& range,
251  const bool debug
252 )
253 {
254  stringList selections(range.size());
255 
256  label nElem = 0;
257 
258  for (int elemI = range.start(); elemI < range.end(); ++elemI)
259  {
260  if (select->GetArraySetting(elemI))
261  {
262  selections[nElem++] = select->GetArrayName(elemI);
263  }
264  }
265 
266  selections.setSize(nElem);
267 
268  // Debugging ...
269  DebugInFunction<< "available(";
270  for (int elemI = range.start(); elemI < range.end(); ++elemI)
271  {
272  DebugInfo<< " \"" << select->GetArrayName(elemI) << "\"";
273  }
274  DebugInfo<< " )\nselected(";
275  forAll(selections, elemI)
276  {
277  DebugInfo<< " " << selections[elemI];
278  }
279  DebugInfo<< " )" << endl;
280 
281  return selections;
282 }
283 
284 
285 void Foam::vtkPVFoam::setSelectedArrayEntries
286 (
287  vtkDataArraySelection* select,
288  const stringList& selections
289 )
290 {
291  const int nElem = select->GetNumberOfArrays();
292  select->DisableAllArrays();
293 
294  // Loop through entries, setting values from selectedEntries
295  for (int elemI=0; elemI < nElem; ++elemI)
296  {
297  string arrayName(select->GetArrayName(elemI));
298 
299  forAll(selections, elemI)
300  {
301  if (selections[elemI] == arrayName)
302  {
303  select->EnableArray(arrayName.c_str());
304  break;
305  }
306  }
307  }
308 }
309 
310 
311 // ************************************************************************* //
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
A class for handling words, derived from string.
Definition: word.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
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
List< string > stringList
A List of strings.
Definition: stringList.H:50