vtkPVFoamFields.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 \*---------------------------------------------------------------------------*/
25 
26 // VTK includes
27 #include "vtkDataArraySelection.h"
28 #include "vtkPolyData.h"
29 #include "vtkUnstructuredGrid.h"
30 
31 // OpenFOAM includes
32 #include "IOobjectList.H"
33 #include "vtkPVFoam.H"
34 #include "vtkPVFoamReader.h"
35 #include "vtkPVFoamVolFields.H"
36 #include "vtkPVFoamPointFields.H"
38 
39 
40 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
41 
42 Foam::IOobjectList Foam::vtkPVFoam::getObjects
43 (
44  const wordHashSet& selected,
45  const fvMesh& mesh,
46  const fileName& instance,
47  const fileName& local
48 )
49 {
50  // If nothing is selected then return an empty list
51  if (selected.empty())
52  {
53  return IOobjectList(0);
54  }
55 
56  // Create the list of objects at the instance
57  IOobjectList objects(mesh, instance, local);
58 
59  // Add any objects from constant that are not already present
60  IOobjectList objectsConstant(mesh, dbPtr_().constant(), local);
61  forAllIter(IOobjectList, objectsConstant, iter)
62  {
63  if (!objects.found(iter.key()))
64  {
65  objects.add
66  (
67  *objectsConstant.HashPtrTable<IOobject>::remove(iter)
68  );
69  }
70  }
71 
72  // Remove everything that is not selected
73  forAllIter(IOobjectList, objects, iter)
74  {
75  if (!selected.found(iter()->name()))
76  {
77  objects.erase(iter);
78  }
79  }
80 
81  return objects;
82 }
83 
84 
85 void Foam::vtkPVFoam::convertFields(vtkMultiBlockDataSet* output)
86 {
87  const fvMesh& mesh = *meshPtr_;
88 
89  wordHashSet selectedFields = getSelected
90  (
91  reader_->GetFieldSelection()
92  );
93 
94  if (selectedFields.empty())
95  {
96  return;
97  }
98 
99  // Get objects (fields) for this time - only keep selected fields
100  // the region name is already in the mesh db
101  IOobjectList objects
102  (
103  getObjects
104  (
105  selectedFields,
106  mesh,
107  dbPtr_().timeName()
108  )
109  );
110 
111  if (objects.empty())
112  {
113  return;
114  }
115 
116  if (debug)
117  {
118  InfoInFunction<< nl
119  << " converting OpenFOAM volume fields" << endl;
120  forAllConstIter(IOobjectList, objects, iter)
121  {
122  Info<< " " << iter()->name()
123  << " == " << iter()->objectPath() << nl;
124  }
125  printMemory();
126  }
127 
128 
129  PtrList<PrimitivePatchInterpolation<primitivePatch>>
130  ppInterpList(mesh.boundaryMesh().size());
131 
132  forAll(ppInterpList, i)
133  {
134  ppInterpList.set
135  (
136  i,
137  new PrimitivePatchInterpolation<primitivePatch>
138  (
139  mesh.boundaryMesh()[i]
140  )
141  );
142  }
143 
144 
145  bool interpFields = reader_->GetInterpolateVolFields();
146 
147  convertVolFields<scalar>
148  (
149  mesh, ppInterpList, objects, interpFields, output
150  );
151  convertVolFields<vector>
152  (
153  mesh, ppInterpList, objects, interpFields, output
154  );
155  convertVolFields<sphericalTensor>
156  (
157  mesh, ppInterpList, objects, interpFields, output
158  );
159  convertVolFields<symmTensor>
160  (
161  mesh, ppInterpList, objects, interpFields, output
162  );
163  convertVolFields<tensor>
164  (
165  mesh, ppInterpList, objects, interpFields, output
166  );
167 
168  convertVolInternalFields<scalar>
169  (
170  mesh, objects, output
171  );
172  convertVolInternalFields<vector>
173  (
174  mesh, objects, output
175  );
176  convertVolInternalFields<sphericalTensor>
177  (
178  mesh, objects, output
179  );
180  convertVolInternalFields<symmTensor>
181  (
182  mesh, objects, output
183  );
184  convertVolInternalFields<tensor>
185  (
186  mesh, objects, output
187  );
188 
189  convertSurfaceFields<scalar>(mesh, objects, output);
190  convertSurfaceFields<vector>(mesh, objects, output);
191  convertSurfaceFields<sphericalTensor>(mesh, objects, output);
192  convertSurfaceFields<symmTensor>(mesh, objects, output);
193  convertSurfaceFields<tensor>(mesh, objects, output);
194 
195  // Construct interpolation on the raw mesh
196  const pointMesh& pMesh = pointMesh::New(mesh);
197 
198  convertPointFields<scalar>(mesh, pMesh, objects, output);
199  convertPointFields<vector>(mesh, pMesh, objects, output);
200  convertPointFields<sphericalTensor>(mesh, pMesh, objects, output);
201  convertPointFields<symmTensor>(mesh, pMesh, objects, output);
202  convertPointFields<tensor>(mesh, pMesh, objects, output);
203 
204  if (debug)
205  {
206  printMemory();
207  }
208 }
209 
210 
211 void Foam::vtkPVFoam::convertLagrangianFields(vtkMultiBlockDataSet* output)
212 {
213  arrayRange& range = arrayRangeLagrangian_;
214  const fvMesh& mesh = *meshPtr_;
215 
216  wordHashSet selectedFields = getSelected
217  (
218  reader_->GetLagrangianFieldSelection()
219  );
220 
221  if (selectedFields.empty())
222  {
223  return;
224  }
225 
226  if (debug)
227  {
228  InfoInFunction << endl;
229  printMemory();
230  }
231 
232  for (int partId = range.start(); partId < range.end(); ++partId)
233  {
234  const word cloudName = getPartName(partId);
235  const label datasetNo = partDataset_[partId];
236 
237  if (!partStatus_[partId] || datasetNo < 0)
238  {
239  continue;
240  }
241 
242 
243  // Get the Lagrangian fields for this time and this cloud
244  // but only keep selected fields
245  // the region name is already in the mesh db
246  IOobjectList objects
247  (
248  getObjects
249  (
250  selectedFields,
251  mesh,
252  dbPtr_().timeName(),
253  cloud::prefix/cloudName
254  )
255  );
256 
257  if (objects.empty())
258  {
259  continue;
260  }
261 
262  if (debug)
263  {
265  << "converting OpenFOAM lagrangian fields" << nl << " ";
266  forAllConstIter(IOobjectList, objects, iter)
267  {
268  Info<< " " << iter()->name()
269  << " == " << iter()->objectPath() << nl;
270  }
271  }
272 
273  convertLagrangianFields<label>
274  (
275  objects, output, datasetNo
276  );
277  convertLagrangianFields<scalar>
278  (
279  objects, output, datasetNo
280  );
281  convertLagrangianFields<vector>
282  (
283  objects, output, datasetNo
284  );
285  convertLagrangianFields<sphericalTensor>
286  (
287  objects, output, datasetNo
288  );
289  convertLagrangianFields<symmTensor>
290  (
291  objects, output, datasetNo
292  );
293  convertLagrangianFields<tensor>
294  (
295  objects, output, datasetNo
296  );
297  }
298 
299  if (debug)
300  {
301  printMemory();
302  }
303 }
304 
305 
306 // ************************************************************************* //
#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
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
static void printMemory()
Simple memory used debugging information.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:208
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
static const char nl
Definition: Ostream.H:260
objects
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:62
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
messageStream Info
#define InfoInFunction
Report an information message using Foam::Info.