1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
27 \*---------------------------------------------------------------------------*/
29 #include "vtkPV3Foam.H"
32 #include "IOobjectList.H"
33 #include "vtkPV3FoamReader.h"
36 #include "vtkDataArraySelection.h"
37 #include "vtkPolyData.h"
38 #include "vtkUnstructuredGrid.h"
40 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42 #include "vtkPV3FoamVolFields.H"
43 #include "vtkPV3FoamPointFields.H"
44 #include "vtkPV3FoamLagrangianFields.H"
47 void Foam::vtkPV3Foam::pruneObjectList
49 IOobjectList& objects,
50 const wordHashSet& selected
53 // hash all the selected field names
59 // only keep selected fields
60 forAllIter(IOobjectList, objects, iter)
62 if (!selected.found(iter()->name()))
70 void Foam::vtkPV3Foam::convertVolFields
72 vtkMultiBlockDataSet* output
75 const fvMesh& mesh = *meshPtr_;
77 wordHashSet selectedFields = getSelected
79 reader_->GetVolFieldSelection()
82 if (!selectedFields.size())
87 // Get objects (fields) for this time - only keep selected fields
88 // the region name is already in the mesh db
89 IOobjectList objects(mesh, dbPtr_().timeName());
90 pruneObjectList(objects, selectedFields);
99 Info<< "<beg> Foam::vtkPV3Foam::convertVolFields" << nl
100 << "converting Foam volume fields" << endl;
101 forAllConstIter(IOobjectList, objects, iter)
103 Info<< " " << iter()->name()
104 << " == " << iter()->objectPath() << nl;
110 PtrList<PrimitivePatchInterpolation<primitivePatch> >
111 ppInterpList(mesh.boundaryMesh().size());
113 forAll(ppInterpList, i)
118 new PrimitivePatchInterpolation<primitivePatch>
120 mesh.boundaryMesh()[i]
126 convertVolFields<scalar>
128 mesh, ppInterpList, objects, output
130 convertVolFields<vector>
132 mesh, ppInterpList, objects, output
134 convertVolFields<sphericalTensor>
136 mesh, ppInterpList, objects, output
138 convertVolFields<symmTensor>
140 mesh, ppInterpList, objects, output
142 convertVolFields<tensor>
144 mesh, ppInterpList, objects, output
149 Info<< "<end> Foam::vtkPV3Foam::convertVolFields" << endl;
155 void Foam::vtkPV3Foam::convertPointFields
157 vtkMultiBlockDataSet* output
160 const fvMesh& mesh = *meshPtr_;
162 wordHashSet selectedFields = getSelected
164 reader_->GetPointFieldSelection()
167 if (!selectedFields.size())
172 // Get objects (fields) for this time - only keep selected fields
173 // the region name is already in the mesh db
174 IOobjectList objects(mesh, dbPtr_().timeName());
175 pruneObjectList(objects, selectedFields);
184 Info<< "<beg> Foam::vtkPV3Foam::convertPointFields" << nl
185 << "converting Foam volume fields" << endl;
186 forAllConstIter(IOobjectList, objects, iter)
188 Info<< " " << iter()->name()
189 << " == " << iter()->objectPath() << nl;
194 // Construct interpolation on the raw mesh
196 // HJ, bug fix? Point mesh handled by objectRegistry
198 const pointMesh& pMesh = pointMesh::New(mesh);
199 // pointMesh pMesh(mesh);
202 convertPointFields<scalar>
204 mesh, pMesh, objects, output
206 convertPointFields<vector>
208 mesh, pMesh, objects, output
210 convertPointFields<sphericalTensor>
212 mesh, pMesh, objects, output
214 convertPointFields<symmTensor>
216 mesh, pMesh, objects, output
218 convertPointFields<tensor>
220 mesh, pMesh, objects, output
225 Info<< "<end> Foam::vtkPV3Foam::convertPointFields" << endl;
231 void Foam::vtkPV3Foam::convertLagrangianFields
233 vtkMultiBlockDataSet* output
236 partInfo& selector = partInfoLagrangian_;
237 const fvMesh& mesh = *meshPtr_;
239 wordHashSet selectedFields = getSelected
241 reader_->GetLagrangianFieldSelection()
244 if (!selectedFields.size())
251 Info<< "<beg> Foam::vtkPV3Foam::convertLagrangianFields" << endl;
255 for (int partId = selector.start(); partId < selector.end(); ++partId)
257 const word cloudName = getPartName(partId);
258 const label datasetNo = partDataset_[partId];
260 if (!partStatus_[partId] || datasetNo < 0)
266 // Get the Lagrangian fields for this time and this cloud
267 // but only keep selected fields
268 // the region name is already in the mesh db
273 cloud::prefix/cloudName
275 pruneObjectList(objects, selectedFields);
284 Info<< "converting Foam lagrangian fields" << nl;
285 forAllConstIter(IOobjectList, objects, iter)
287 Info<< " " << iter()->name()
288 << " == " << iter()->objectPath() << nl;
292 convertLagrangianFields<label>
294 objects, output, datasetNo
296 convertLagrangianFields<scalar>
298 objects, output, datasetNo
300 convertLagrangianFields<vector>
302 objects, output, datasetNo
304 convertLagrangianFields<sphericalTensor>
306 objects, output, datasetNo
308 convertLagrangianFields<symmTensor>
310 objects, output, datasetNo
312 convertLagrangianFields<tensor>
314 objects, output, datasetNo
320 Info<< "<end> Foam::vtkPV3Foam::convertLagrangianFields" << endl;
326 // ************************************************************************* //