patchWriter.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011 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 "writeFuns.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  patchIDs_(patchIDs),
45  // Use binary mode in case we write binary.
46  // Causes windows reading to fail if we don't
47  os_(fName.c_str(),
48  ios_base::out|ios_base::binary)
49 {
50  const fvMesh& mesh = vMesh_.mesh();
52 
53  // Write header
54  if (patchIDs_.size() == 1)
55  {
56  writeFuns::writeHeader(os_, binary_, patches[patchIDs_[0]].name());
57  }
58  else
59  {
60  writeFuns::writeHeader(os_, binary_, "patches");
61  }
62  os_ << "DATASET POLYDATA" << std::endl;
63 
64  // Write topology
65  nPoints_ = 0;
66  nFaces_ = 0;
67  label nFaceVerts = 0;
68 
69  forAll(patchIDs_, i)
70  {
71  const polyPatch& pp = patches[patchIDs_[i]];
72 
73  nPoints_ += pp.nPoints();
74  nFaces_ += pp.size();
75 
76  forAll(pp, faceI)
77  {
78  nFaceVerts += pp[faceI].size() + 1;
79  }
80  }
81 
82  os_ << "POINTS " << nPoints_ << " double" << std::endl;
83 
84  DynamicList<doubleScalar> ptField(3*nPoints_);
85 
86  forAll(patchIDs_, i)
87  {
88  const polyPatch& pp = patches[patchIDs_[i]];
89 
90  writeFuns::insert(pp.localPoints(), ptField);
91  }
92  writeFuns::write(os_, binary_, ptField);
93 
94  os_ << "POLYGONS " << nFaces_ << ' ' << nFaceVerts << std::endl;
95 
96  DynamicList<label> vertLabels(nFaceVerts);
97 
98  label offset = 0;
99 
100  forAll(patchIDs_, i)
101  {
102  const polyPatch& pp = patches[patchIDs_[i]];
103 
104  forAll(pp, faceI)
105  {
106  const face& f = pp.localFaces()[faceI];
107 
108  vertLabels.append(f.size());
109  writeFuns::insert(f + offset, vertLabels);
110  }
111  offset += pp.nPoints();
112  }
113  writeFuns::write(os_, binary_, vertLabels);
114 }
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 {
121  const fvMesh& mesh = vMesh_.mesh();
122 
124 
125  os_ << "patchID 1 " << nFaces_ << " double" << std::endl;
126 
127  forAll(patchIDs_, i)
128  {
129  label patchI = patchIDs_[i];
130 
131  const polyPatch& pp = mesh.boundaryMesh()[patchI];
132 
133  if (!isA<emptyPolyPatch>(pp))
134  {
135  writeFuns::insert(scalarField(pp.size(), patchI), fField);
136  }
137  }
138  writeFuns::write(os_, binary_, fField);
139 }
140 
141 
142 // ************************************************************************* //
patchWriter.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
writeFuns.H
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::vtkMesh
Encapsulation of VTK mesh data. Holds mesh or meshsubset and polyhedral-cell decomposition on it.
Definition: vtkMesh.H:52
Foam::patchWriter::writePatchIDs
void writePatchIDs()
Write cellIDs.
Definition: patchWriter.C:119
Foam::label
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
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::patchWriter::nFaces_
label nFaces_
Definition: patchWriter.H:75
Foam::writeFuns::write
static void write(std::ostream &, const bool, DynamicList< floatScalar > &)
Write floats ascii or binary.
Definition: writeFuns.C:107
Foam::patchWriter::patchIDs_
const labelList patchIDs_
Definition: patchWriter.H:69
Foam::PrimitivePatch::nPoints
label nPoints() const
Return number of points supporting patch faces.
Definition: PrimitivePatchTemplate.H:293
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::vtkMesh::mesh
const fvMesh & mesh() const
Access either mesh or submesh.
Definition: vtkMesh.H:119
f
labelList f(nPoints)
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::patchWriter::vMesh_
const vtkMesh & vMesh_
Definition: patchWriter.H:61
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatchTemplate.C:372
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::writeFuns::insert
static void insert(const point &, DynamicList< floatScalar > &dest)
Append point to given DynamicList.
Definition: writeFuns.C:170
write
Tcoeff write()
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
insert
timeIndices insert(timeIndex, timeDirs[timeI].value())
Foam::patchWriter::binary_
const bool binary_
Definition: patchWriter.H:63
Foam::patchWriter::patchWriter
patchWriter(const vtkMesh &, const bool binary, const bool nearCellValue, const fileName &, const labelList &patchIDs)
Construct from components.
Definition: patchWriter.C:32
Foam::patchWriter::os_
std::ofstream os_
Definition: patchWriter.H:71