indexedOctreeBase.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) 2011-2012 OpenFOAM Foundation
9  Copyright (C) 2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "indexedOctree.H"
30 #include "treeBoundBox.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(indexedOctreeBase, 0);
37 
38  scalar indexedOctreeBase::perturbTol_ = 10*SMALL;
39 }
40 
41 
42 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
43 
45 (
46  Ostream& os,
47  const treeBoundBox& bb,
48  label& vertIndex,
49  const bool writeLinesOnly
50 )
51 {
52  // Annotate with '#box' which can be grep'd for later
53  os << "#box" << nl;
54 
55  pointField pts(bb.points());
56 
57  for (const point& p : pts)
58  {
59  os << "v " << p.x() << ' ' << p.y() << ' ' << p.z() << nl;
60  }
61 
62  if (writeLinesOnly)
63  {
64  for (const edge& e : treeBoundBox::edges)
65  {
66  os << "l " << e[0] + vertIndex + 1
67  << ' ' << e[1] + vertIndex + 1 << nl;
68  }
69  }
70  else
71  {
72  for (const face& f : treeBoundBox::faces)
73  {
74  os << 'f';
75  for (const label fpi : f)
76  {
77  os << ' ' << fpi + vertIndex + 1;
78  }
79  os << nl;
80  }
81  }
82 
83  vertIndex += pts.size();
84 }
85 
86 
87 // ************************************************************************* //
static void writeOBJ(Ostream &os, const treeBoundBox &bb, label &vertIndex, const bool writeLinesOnly=false)
Write treeBoundBox in OBJ format.
static scalar perturbTol_
Relative perturbation tolerance.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
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
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)
defineTypeNameAndDebug(combustionModel, 0)
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
Namespace for OpenFOAM.
const pointField & pts