initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / graphics / PV3FoamReader / vtkPV3Foam / vtkPV3FoamMeshLagrangian.C
blob07d322d54250876e929526d46274ecf8e3bece92
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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
25 Description
27 \*---------------------------------------------------------------------------*/
29 #include "vtkPV3Foam.H"
31 // Foam includes
32 #include "Cloud.H"
33 #include "fvMesh.H"
34 #include "IOobjectList.H"
35 #include "passiveParticle.H"
36 #include "vtkPV3FoamPoints.H"
38 // VTK includes
39 #include "vtkCellArray.h"
40 #include "vtkPoints.h"
41 #include "vtkPolyData.h"
43 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
45 vtkPolyData* Foam::vtkPV3Foam::lagrangianVTKMesh
47     const fvMesh& mesh,
48     const word& cloudName
51     vtkPolyData* vtkmesh = NULL;
53     if (debug)
54     {
55         Info<< "<beg> Foam::vtkPV3Foam::lagrangianVTKMesh - timePath "
56             << mesh.time().timePath()/cloud::prefix/cloudName << endl;
57         printMemory();
58     }
61     // the region name is already in the mesh db
62     IOobjectList sprayObjs
63     (
64         mesh,
65         mesh.time().timeName(),
66         cloud::prefix/cloudName
67     );
69     IOobject* positionsPtr = sprayObjs.lookup("positions");
70     if (positionsPtr)
71     {
72         Cloud<passiveParticle> parcels(mesh, cloudName, false);
74         if (debug)
75         {
76             Info<< "cloud with " << parcels.size() << " parcels" << endl;
77         }
79         vtkmesh = vtkPolyData::New();
80         vtkPoints* vtkpoints = vtkPoints::New();
81         vtkCellArray* vtkcells = vtkCellArray::New();
83         vtkpoints->Allocate( parcels.size() );
84         vtkcells->Allocate( parcels.size() );
86         vtkIdType particleId = 0;
87         forAllConstIter(Cloud<passiveParticle>, parcels, iter)
88         {
89             vtkPV3FoamInsertNextPoint(vtkpoints, iter().position());
91             vtkcells->InsertNextCell(1, &particleId);
92             particleId++;
93         }
95         vtkmesh->SetPoints(vtkpoints);
96         vtkpoints->Delete();
98         vtkmesh->SetVerts(vtkcells);
99         vtkcells->Delete();
100     }
102     if (debug)
103     {
104         Info<< "<end> Foam::vtkPV3Foam::lagrangianVTKMesh" << endl;
105         printMemory();
106     }
108     return vtkmesh;
112 // ************************************************************************* //