vtkPV4FoamFields.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-2014 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 "vtkPV4Foam.H"
27 
28 // OpenFOAM includes
29 #include "IOobjectList.H"
30 #include "vtkPV4FoamReader.h"
31 
32 // VTK includes
33 #include "vtkDataArraySelection.h"
34 #include "vtkPolyData.h"
35 #include "vtkUnstructuredGrid.h"
36 
37 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38 
39 #include "vtkPV4FoamVolFields.H"
40 #include "vtkPV4FoamPointFields.H"
42 
43 
45 (
46  IOobjectList& objects,
47  const wordHashSet& selected
48 )
49 {
50  // hash all the selected field names
51  if (selected.empty())
52  {
53  objects.clear();
54  }
55 
56  // only keep selected fields
57  forAllIter(IOobjectList, objects, iter)
58  {
59  if (!selected.found(iter()->name()))
60  {
61  objects.erase(iter);
62  }
63  }
64 }
65 
66 
68 (
69  vtkMultiBlockDataSet* output
70 )
71 {
72  const fvMesh& mesh = *meshPtr_;
73 
74  wordHashSet selectedFields = getSelected
75  (
76  reader_->GetVolFieldSelection()
77  );
78 
79  if (selectedFields.empty())
80  {
81  return;
82  }
83 
84  // Get objects (fields) for this time - only keep selected fields
85  // the region name is already in the mesh db
86  IOobjectList objects(mesh, dbPtr_().timeName());
87  pruneObjectList(objects, selectedFields);
88 
89  if (objects.empty())
90  {
91  return;
92  }
93 
94  if (debug)
95  {
96  Info<< "<beg> Foam::vtkPV4Foam::convertVolFields" << nl
97  << "converting OpenFOAM volume fields" << endl;
98  forAllConstIter(IOobjectList, objects, iter)
99  {
100  Info<< " " << iter()->name()
101  << " == " << iter()->objectPath() << nl;
102  }
103  printMemory();
104  }
105 
106 
107  PtrList<PrimitivePatchInterpolation<primitivePatch> >
108  ppInterpList(mesh.boundaryMesh().size());
109 
110  forAll(ppInterpList, i)
111  {
112  ppInterpList.set
113  (
114  i,
115  new PrimitivePatchInterpolation<primitivePatch>
116  (
117  mesh.boundaryMesh()[i]
118  )
119  );
120  }
121 
122 
123  bool interpFields = reader_->GetInterpolateVolFields();
124 
125  convertVolFields<scalar>
126  (
127  mesh, ppInterpList, objects, interpFields, output
128  );
129  convertVolFields<vector>
130  (
131  mesh, ppInterpList, objects, interpFields, output
132  );
133  convertVolFields<sphericalTensor>
134  (
135  mesh, ppInterpList, objects, interpFields, output
136  );
137  convertVolFields<symmTensor>
138  (
139  mesh, ppInterpList, objects, interpFields, output
140  );
141  convertVolFields<tensor>
142  (
143  mesh, ppInterpList, objects, interpFields, output
144  );
145 
146  convertDimFields<scalar>
147  (
148  mesh, ppInterpList, objects, interpFields, output
149  );
150  convertDimFields<vector>
151  (
152  mesh, ppInterpList, objects, interpFields, output
153  );
154  convertDimFields<sphericalTensor>
155  (
156  mesh, ppInterpList, objects, interpFields, output
157  );
158  convertDimFields<symmTensor>
159  (
160  mesh, ppInterpList, objects, interpFields, output
161  );
162  convertDimFields<tensor>
163  (
164  mesh, ppInterpList, objects, interpFields, output
165  );
166 
167  if (debug)
168  {
169  Info<< "<end> Foam::vtkPV4Foam::convertVolFields" << endl;
170  printMemory();
171  }
172 }
173 
174 
176 (
177  vtkMultiBlockDataSet* output
178 )
179 {
180  const fvMesh& mesh = *meshPtr_;
181 
182  wordHashSet selectedFields = getSelected
183  (
184  reader_->GetPointFieldSelection()
185  );
186 
187  if (selectedFields.empty())
188  {
189  if (debug)
190  {
191  Info<< "no point fields selected" << endl;
192  }
193  return;
194  }
195 
196  // Get objects (fields) for this time - only keep selected fields
197  // the region name is already in the mesh db
198  IOobjectList objects(mesh, dbPtr_().timeName());
199  pruneObjectList(objects, selectedFields);
200 
201  if (objects.empty())
202  {
203  return;
204  }
205 
206  if (debug)
207  {
208  Info<< "<beg> Foam::vtkPV4Foam::convertPointFields" << nl
209  << "converting OpenFOAM volume fields -> point fields" << endl;
210  forAllConstIter(IOobjectList, objects, iter)
211  {
212  Info<< " " << iter()->name()
213  << " == " << iter()->objectPath() << nl;
214  }
215  printMemory();
216  }
217 
218  // Construct interpolation on the raw mesh
219  const pointMesh& pMesh = pointMesh::New(mesh);
220 
221 
222  convertPointFields<scalar>
223  (
224  mesh, pMesh, objects, output
225  );
226  convertPointFields<vector>
227  (
228  mesh, pMesh, objects, output
229  );
230  convertPointFields<sphericalTensor>
231  (
232  mesh, pMesh, objects, output
233  );
234  convertPointFields<symmTensor>
235  (
236  mesh, pMesh, objects, output
237  );
238  convertPointFields<tensor>
239  (
240  mesh, pMesh, objects, output
241  );
242 
243  if (debug)
244  {
245  Info<< "<end> Foam::vtkPV4Foam::convertPointFields" << endl;
246  printMemory();
247  }
248 }
249 
250 
252 (
253  vtkMultiBlockDataSet* output
254 )
255 {
256  arrayRange& range = arrayRangeLagrangian_;
257  const fvMesh& mesh = *meshPtr_;
258 
259  wordHashSet selectedFields = getSelected
260  (
261  reader_->GetLagrangianFieldSelection()
262  );
263 
264  if (selectedFields.empty())
265  {
266  return;
267  }
268 
269  if (debug)
270  {
271  Info<< "<beg> Foam::vtkPV4Foam::convertLagrangianFields" << endl;
272  printMemory();
273  }
274 
275  for (int partId = range.start(); partId < range.end(); ++partId)
276  {
277  const word cloudName = getPartName(partId);
278  const label datasetNo = partDataset_[partId];
279 
280  if (!partStatus_[partId] || datasetNo < 0)
281  {
282  continue;
283  }
284 
285 
286  // Get the Lagrangian fields for this time and this cloud
287  // but only keep selected fields
288  // the region name is already in the mesh db
289  IOobjectList objects
290  (
291  mesh,
292  dbPtr_().timeName(),
294  );
295  pruneObjectList(objects, selectedFields);
296 
297  if (objects.empty())
298  {
299  continue;
300  }
301 
302  if (debug)
303  {
304  Info<< "converting OpenFOAM lagrangian fields" << nl;
305  forAllConstIter(IOobjectList, objects, iter)
306  {
307  Info<< " " << iter()->name()
308  << " == " << iter()->objectPath() << nl;
309  }
310  }
311 
312  convertLagrangianFields<label>
313  (
314  objects, output, datasetNo
315  );
316  convertLagrangianFields<scalar>
317  (
318  objects, output, datasetNo
319  );
320  convertLagrangianFields<vector>
321  (
322  objects, output, datasetNo
323  );
324  convertLagrangianFields<sphericalTensor>
325  (
326  objects, output, datasetNo
327  );
328  convertLagrangianFields<symmTensor>
329  (
330  objects, output, datasetNo
331  );
332  convertLagrangianFields<tensor>
333  (
334  objects, output, datasetNo
335  );
336  }
337 
338  if (debug)
339  {
340  Info<< "<end> Foam::vtkPV4Foam::convertLagrangianFields" << endl;
341  printMemory();
342  }
343 }
344 
345 
346 // ************************************************************************* //
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
vtkPV4Foam.H
IOobjectList.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
vtkPV4FoamLagrangianFields.H
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
vtkPV4FoamVolFields.H
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::vtkPV3Readers::getSelected
wordHashSet getSelected(vtkDataArraySelection *select)
Retrieve the current selections as a wordHashSet.
Definition: vtkPV3Readers.C:187
Foam::label
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
Foam::vtkPV4Foam::convertLagrangianFields
void convertLagrangianFields(vtkMultiBlockDataSet *)
Convert Lagrangian fields.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
vtkPV4FoamPointFields.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::cloud::prefix
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
range
scalar range
Definition: LISASMDCalcMethod1.H:12
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
timeName
word timeName
Definition: getTimeIndex.H:3
cloudName
const word cloudName(propsDict.lookup("cloudName"))
Foam::wordHashSet
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:207
Foam::vtkPV4Foam::pruneObjectList
static void pruneObjectList(IOobjectList &, const wordHashSet &)
Only keep what is listed in hashSet.
Foam::vtkPV4Foam::convertPointFields
void convertPointFields(vtkMultiBlockDataSet *)
Convert point fields.
Foam::vtkPV4Foam::convertVolFields
void convertVolFields(vtkMultiBlockDataSet *)
Convert volume fields.
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47