DOC: Corrected class names in the file descriptions
[freefoam.git] / applications / utilities / postProcessing / dataConversion / foamToVTK / lagrangianWriter.C
blob435379c498eacf1cd34d1faced108d7351858bc9
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "lagrangianWriter.H"
27 #include "writeFuns.H"
28 #include <lagrangian/Cloud.H>
29 #include <lagrangian/passiveParticle.H>
31 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
33 // Construct from components
34 Foam::lagrangianWriter::lagrangianWriter
36     const vtkMesh& vMesh,
37     const bool binary,
38     const fileName& fName,
39     const word& cloudName,
40     const bool dummyCloud
43     vMesh_(vMesh),
44     binary_(binary),
45     fName_(fName),
46     cloudName_(cloudName),
47     os_(fName.c_str())
49     const fvMesh& mesh = vMesh_.mesh();
51     // Write header
52     writeFuns::writeHeader(os_, binary_, mesh.time().caseName());
53     os_ << "DATASET POLYDATA" << std::endl;
55     if (dummyCloud)
56     {
57         nParcels_ = 0;
59         os_ << "POINTS " << nParcels_ << " float" << std::endl;
60     }
61     else
62     {
63         Cloud<passiveParticle> parcels(mesh, cloudName_, false);
65         nParcels_ = parcels.size();
67         os_ << "POINTS " << nParcels_ << " float" << std::endl;
69         DynamicList<floatScalar> partField(3*parcels.size());
71         forAllConstIter(Cloud<passiveParticle>, parcels, elmnt)
72         {
73             writeFuns::insert(elmnt().position(), partField);
74         }
75         writeFuns::write(os_, binary_, partField);
76     }
80 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
82 void Foam::lagrangianWriter::writeParcelHeader(const label nFields)
84     os_ << "POINT_DATA " << nParcels_ << std::endl
85         << "FIELD attributes " << nFields
86         << std::endl;
90 // ************************ vim: set sw=4 sts=4 et: ************************ //