vtkPV3Readers.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 "vtkPV3Readers.H"
30 
31 // OpenFOAM includes
32 #include "IFstream.H"
33 
34 // VTK includes
35 #include "vtkDataArraySelection.h"
36 #include "vtkDataSet.h"
37 #include "vtkMultiBlockDataSet.h"
38 #include "vtkInformation.h"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 defineTypeNameAndDebug(vtkPV3Readers, 0);
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
53  // Extract up to the first non-word characters
54  inline word getFirstWord(const char* str)
55  {
56  if (str)
57  {
58  label n = 0;
59  while (str[n] && word::valid(str[n]))
60  {
61  ++n;
62  }
63  return word(str, n, true);
64  }
65  else
66  {
67  return word::null;
68  }
69 
70  }
72 
73 } // End namespace Foam
74 
75 
76 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
77 
79 (
80  vtkMultiBlockDataSet* output,
81  vtkDataSet* dataset,
82  const partInfo& selector,
83  const label datasetNo,
84  const std::string& datasetName
85 )
86 {
87  const int blockNo = selector.block();
88 
89  vtkDataObject* blockDO = output->GetBlock(blockNo);
90  vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
91 
92  if (!block)
93  {
94  if (blockDO)
95  {
97  << "Block already has a vtkDataSet assigned to it"
98  << endl;
99  return;
100  }
101 
102  block = vtkMultiBlockDataSet::New();
103  output->SetBlock(blockNo, block);
104  block->Delete();
105  }
106 
107  if (debug)
108  {
109  Info<< "block[" << blockNo << "] has "
110  << block->GetNumberOfBlocks()
111  << " datasets prior to adding set " << datasetNo
112  << " with name: " << datasetName << endl;
113  }
114 
115  block->SetBlock(datasetNo, dataset);
116 
117  // name the block when assigning dataset 0
118  if (datasetNo == 0)
119  {
120  output->GetMetaData(blockNo)->Set
121  (
122  vtkCompositeDataSet::NAME(),
123  selector.name()
124  );
125  }
126 
127  if (datasetName.size())
128  {
129  block->GetMetaData(datasetNo)->Set
130  (
131  vtkCompositeDataSet::NAME(),
132  datasetName.c_str()
133  );
134  }
135 }
136 
137 
139 (
140  vtkMultiBlockDataSet* output,
141  const partInfo& selector,
142  const label datasetNo
143 )
144 {
145  const int blockNo = selector.block();
146 
147  vtkDataObject* blockDO = output->GetBlock(blockNo);
148  vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
149 
150  if (block)
151  {
152  return vtkDataSet::SafeDownCast(block->GetBlock(datasetNo));
153  }
154 
155  return 0;
156 }
157 
158 
159 // ununsed at the moment
161 (
162  vtkMultiBlockDataSet* output,
163  const partInfo& selector
164 )
165 {
166  const int blockNo = selector.block();
167 
168  vtkDataObject* blockDO = output->GetBlock(blockNo);
169  vtkMultiBlockDataSet* block = vtkMultiBlockDataSet::SafeDownCast(blockDO);
170  if (block)
171  {
172  return block->GetNumberOfBlocks();
173  }
174 
175  return 0;
176 }
177 
178 
179 // Foam::word Foam::vtkPV3Readers::getPartName(int partId)
180 // {
181 // return getFirstWord(reader_->GetPartArrayName(partId));
182 // }
183 
184 
186 (
187  vtkDataArraySelection* select
188 )
189 {
190  int nElem = select->GetNumberOfArrays();
191  wordHashSet selections(2*nElem);
192 
193  for (int elemI=0; elemI < nElem; ++elemI)
194  {
195  if (select->GetArraySetting(elemI))
196  {
197  selections.insert(getFirstWord(select->GetArrayName(elemI)));
198  }
199  }
200 
201  return selections;
202 }
203 
204 
206 (
207  vtkDataArraySelection* select,
208  const partInfo& selector
209 )
210 {
211  int nElem = select->GetNumberOfArrays();
212  wordHashSet selections(2*nElem);
213 
214  for (int elemI = selector.start(); elemI < selector.end(); ++elemI)
215  {
216  if (select->GetArraySetting(elemI))
217  {
218  selections.insert(getFirstWord(select->GetArrayName(elemI)));
219  }
220  }
221 
222  return selections;
223 }
224 
225 
227 (
228  vtkDataArraySelection* select
229 )
230 {
231  stringList selections(select->GetNumberOfArrays());
232  label nElem = 0;
233 
234  forAll(selections, elemI)
235  {
236  if (select->GetArraySetting(elemI))
237  {
238  selections[nElem++] = select->GetArrayName(elemI);
239  }
240  }
241  selections.setSize(nElem);
242 
243 
244  if (debug)
245  {
246  label nElem = select->GetNumberOfArrays();
247  Info<< "available(";
248  for (int elemI = 0; elemI < nElem; ++elemI)
249  {
250  Info<< " \"" << select->GetArrayName(elemI) << "\"";
251  }
252  Info<< " )\nselected(";
253 
254  forAll(selections, elemI)
255  {
256  Info<< " " << selections[elemI];
257  }
258  Info<< " )\n";
259  }
260 
261  return selections;
262 }
263 
264 
266 (
267  vtkDataArraySelection* select,
268  const partInfo& selector
269 )
270 {
271  stringList selections(selector.size());
272  label nElem = 0;
273 
274  for (int elemI = selector.start(); elemI < selector.end(); ++elemI)
275  {
276  if (select->GetArraySetting(elemI))
277  {
278  selections[nElem++] = select->GetArrayName(elemI);
279  }
280  }
281  selections.setSize(nElem);
282 
283 
284  if (debug)
285  {
286  Info<< "available(";
287  for (int elemI = selector.start(); elemI < selector.end(); ++elemI)
288  {
289  Info<< " \"" << select->GetArrayName(elemI) << "\"";
290  }
291  Info<< " )\nselected(";
292 
293  forAll(selections, elemI)
294  {
295  Info<< " " << selections[elemI];
296  }
297  Info<< " )\n";
298  }
299 
300  return selections;
301 }
302 
303 
305 (
306  vtkDataArraySelection* select,
307  const stringList& selections
308 )
309 {
310  const int nElem = select->GetNumberOfArrays();
311  select->DisableAllArrays();
312 
313  // Loop through entries, setting values from selectedEntries
314  for (int elemI=0; elemI < nElem; ++elemI)
315  {
316  string arrayName(select->GetArrayName(elemI));
317 
318  forAll(selections, elemI)
319  {
320  if (selections[elemI] == arrayName)
321  {
322  select->EnableArray(arrayName.c_str());
323  break;
324  }
325  }
326  }
327 }
328 
329 
330 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
331 
332 
333 // ************************************************************************* //
void AddToBlock(vtkMultiBlockDataSet *output, vtkDataSet *dataset, const partInfo &selector, const label datasetNo, const std::string &datasetName)
Convenience method use to convert the readers from VTK 5.
Definition: vtkPV3Readers.C:80
static bool valid(char)
Is this character valid for a word.
Definition: wordI.H:117
A HashTable with keys but without contents.
Definition: HashSet.H:59
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
#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:253
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
int block() const
Return the block holding these datasets.
Definition: vtkPV3Readers.H:93
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Bookkeeping for GUI checklists and the multi-block organization.
Definition: vtkPV3Readers.H:75
stringList getSelectedArrayEntries(vtkDataArraySelection *)
Retrieve the current selections.
label GetNumberOfDataSets(vtkMultiBlockDataSet *output, const partInfo &selector)
Convenience method use to convert the readers from VTK 5.
static const word null
An empty word.
Definition: word.H:77
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:63
defineTypeNameAndDebug(combustionModel, 0)
wordHashSet getSelected(vtkDataArraySelection *select)
Retrieve the current selections as a wordHashSet.
void setSize(const label)
Reset size of List.
Definition: List.C:295
vtkDataSet * GetDataSetFromBlock(vtkMultiBlockDataSet *output, const partInfo &selector, const label datasetNo)
Convenience method use to convert the readers from VTK 5.
void setSelectedArrayEntries(vtkDataArraySelection *, const stringList &)
Set selection(s)
messageStream Info
label n
const char * name() const
Namespace for OpenFOAM.