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