patchWriter.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-2023 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 "patchWriter.H"
27 #include "vtkWriteFieldOps.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const vtkMesh& vMesh,
34  const bool binary,
35  const bool nearCellValue,
36  const fileName& fName,
37  const labelList& patchIDs
38 )
39 :
40  vMesh_(vMesh),
41  binary_(binary),
42  nearCellValue_(nearCellValue),
43  fName_(fName),
44  patchIndices_(patchIDs),
45  os_(fName.c_str())
46 {
47  const fvMesh& mesh = vMesh_.mesh();
48  const polyBoundaryMesh& patches = mesh.boundaryMesh();
49 
50  // Write header
51  if (patchIndices_.size() == 1)
52  {
54  (
55  os_,
56  binary_, patches[patchIndices_[0]].name()
57  );
58  }
59  else
60  {
61  vtkWriteOps::writeHeader(os_, binary_, "patches");
62  }
63  os_ << "DATASET POLYDATA" << std::endl;
64 
65  // Write topology
66  nPoints_ = 0;
67  nFaces_ = 0;
68  label nFaceVerts = 0;
69 
70  forAll(patchIndices_, i)
71  {
72  const polyPatch& pp = patches[patchIndices_[i]];
73 
74  nPoints_ += pp.nPoints();
75  nFaces_ += pp.size();
76 
77  forAll(pp, facei)
78  {
79  nFaceVerts += pp[facei].size() + 1;
80  }
81  }
82 
83  os_ << "POINTS " << nPoints_ << " float" << std::endl;
84 
85  DynamicList<floatScalar> ptField(3*nPoints_);
86 
87  forAll(patchIndices_, i)
88  {
89  const polyPatch& pp = patches[patchIndices_[i]];
90 
91  vtkWriteOps::insert(pp.localPoints(), ptField);
92  }
93  vtkWriteOps::write(os_, binary_, ptField);
94 
95  os_ << "POLYGONS " << nFaces_ << ' ' << nFaceVerts << std::endl;
96 
97  DynamicList<label> vertLabels(nFaceVerts);
98 
99  label offset = 0;
100 
101  forAll(patchIndices_, i)
102  {
103  const polyPatch& pp = patches[patchIndices_[i]];
104 
105  forAll(pp, facei)
106  {
107  const face& f = pp.localFaces()[facei];
108 
109  vertLabels.append(f.size());
110  vtkWriteOps::insert(f + offset, vertLabels);
111  }
112  offset += pp.nPoints();
113  }
114  vtkWriteOps::write(os_, binary_, vertLabels);
115 }
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
121 {
122  const fvMesh& mesh = vMesh_.mesh();
123 
124  DynamicList<floatScalar> fField(nFaces_);
125 
126  os_ << "patchID 1 " << nFaces_ << " float" << std::endl;
127 
128  forAll(patchIndices_, i)
129  {
130  label patchi = patchIndices_[i];
131 
132  const polyPatch& pp = mesh.boundaryMesh()[patchi];
133 
134  if (!isA<emptyPolyPatch>(pp))
135  {
136  vtkWriteOps::insert(scalarField(pp.size(), patchi), fField);
137  }
138  }
139  vtkWriteOps::write(os_, binary_, fField);
140 }
141 
142 
143 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
patchWriter(const vtkMesh &, const bool binary, const bool nearCellValue, const fileName &, const labelList &patchIDs)
Construct from components.
void writePatchIndices()
Write cellIDs.
const fvMesh & mesh() const
Access either mesh or submesh.
Definition: vtkMesh.H:122
label patchi
const fvPatchList & patches
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.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void offset(label &lst, const label o)
labelList f(nPoints)