lagrangianWriter.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-2022 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 "lagrangianWriter.H"
27 #include "vtkWriteFieldOps.H"
28 #include "Cloud.H"
29 #include "passiveParticle.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 (
35  const vtkMesh& vMesh,
36  const bool binary,
37  const fileName& fName,
38  const word& cloudName,
39  const bool dummyCloud
40 )
41 :
42  vMesh_(vMesh),
43  binary_(binary),
44  fName_(fName),
45  cloudName_(cloudName),
46  os_(fName.c_str())
47 {
48  const fvMesh& mesh = vMesh_.mesh();
49 
50  // Write header
51  vtkWriteOps::writeHeader(os_, binary_, mesh.time().caseName());
52  os_ << "DATASET POLYDATA" << std::endl;
53 
54  if (dummyCloud)
55  {
56  nParcels_ = 0;
57 
58  os_ << "POINTS " << nParcels_ << " float" << std::endl;
59 
60  os_ << "VERTICES " << nParcels_ << ' ' << 2*nParcels_ << std::endl;
61  }
62  else
63  {
64  Cloud<passiveParticle> parcels(mesh, cloudName_, false);
65 
66  nParcels_ = parcels.size();
67 
68  os_ << "POINTS " << nParcels_ << " float" << std::endl;
69 
70  DynamicList<floatScalar> partField(3*parcels.size());
71  forAllConstIter(Cloud<passiveParticle>, parcels, elmnt)
72  {
73  vtkWriteOps::insert(elmnt().position(mesh), partField);
74  }
75  vtkWriteOps::write(os_, binary_, partField);
76 
77  os_ << "VERTICES " << nParcels_ << ' ' << 2*nParcels_ << std::endl;
78 
79  DynamicList<label> vertexPoints(2*parcels.size());
80  forAll(parcels, parceli)
81  {
82  vertexPoints.append(1);
83  vertexPoints.append(parceli);
84  }
85  vtkWriteOps::write(os_, binary, vertexPoints);
86  }
87 }
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  os_ << "POINT_DATA " << nParcels_ << std::endl
95  << "FIELD attributes " << nFields
96  << std::endl;
97 }
98 
99 
100 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
void writeFieldsHeader(const label nFields)
Write the fields header specifying the number of fields.
lagrangianWriter(const vtkMesh &, const bool binary, const fileName &, const word &, const bool dummyCloud)
Construct from components.
void insert(const scalar, DynamicList< floatScalar > &)
Append scalar to given DynamicList.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const word cloudName(propsDict.lookup("cloudName"))