AABBTreeBase.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "AABBTree.H"
29 #include "treeBoundBox.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 Foam::scalar Foam::AABBTreeBase::tolerance_ = 1e-4;
34 
35 
36 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
37 
39 (
40  Ostream& os,
41  const treeBoundBox& bb,
42  label& vertIndex,
43  const bool writeLinesOnly
44 )
45 {
46  // Annotate with '#box' which can be grep'd for later
47  os << "#box" << nl;
48 
49  pointField pts(bb.points());
50 
51  for (const point& p : pts)
52  {
53  os << "v " << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
54  }
55 
56  if (writeLinesOnly)
57  {
58  for (const edge& e : treeBoundBox::edges)
59  {
60  os << "l " << e[0] + vertIndex + 1
61  << ' ' << e[1] + vertIndex + 1 << nl;
62  }
63  }
64  else
65  {
66  for (const face& f : treeBoundBox::faces)
67  {
68  os << 'f';
69  for (const label fpi : f)
70  {
71  os << ' ' << fpi + vertIndex + 1;
72  }
73  os << nl;
74  }
75  }
76 
77  vertIndex += pts.size();
78 }
79 
80 
81 // ************************************************************************* //
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
static void writeOBJ(Ostream &os, const treeBoundBox &bb, label &vertIndex, const bool writeLinesOnly=false)
Write treeBoundBox in OBJ format.
Definition: AABBTreeBase.C:32
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static scalar tolerance_
Relative tolerance.
Definition: AABBTree.H:80
static const faceList faces
Face to point addressing, using octant corner points.
Definition: treeBoundBox.H:174
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
static const edgeList edges
Edge to point addressing, using octant corner points.
Definition: treeBoundBox.H:179
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
labelList f(nPoints)
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:87
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
volScalarField & p
const pointField & pts