PolyhedronReaderTemplates.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) 2015 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 "PolyhedronReader.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class HDS>
32 (
33  const triSurface& s
34 )
35 :
36  s_(s)
37 {}
38 
39 
40 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
41 
42 template<class HDS>
44 {
45  // Postcondition: hds is a valid polyhedral surface.
46  CGAL::Polyhedron_incremental_builder_3<HDS> B(hds, true);
47 
48  B.begin_surface(s_.nPoints(), s_.size());
49 
50  typedef typename HDS::Vertex Vertex;
51  typedef typename Vertex::Point Point;
52 
53  const Foam::pointField& pts = s_.points();
54  forAll(pts, i)
55  {
56  const Foam::point& pt = pts[i];
57  B.add_vertex(Point(pt.x(), pt.y(), pt.z()));
58  }
59  forAll(s_, i)
60  {
61  const Foam::labelledTri& t = s_[i];
62  B.begin_facet();
63  B.add_vertex_to_facet(t[0]);
64  B.add_vertex_to_facet(t[1]);
65  B.add_vertex_to_facet(t[2]);
66  B.end_facet();
67  }
68  B.end_surface();
69 }
70 
71 
72 // ************************************************************************* //
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Vertex
graph_traits< Graph >::vertex_descriptor Vertex
Definition: SloanRenumber.C:72
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::PolyhedronReader::Build_triangle::Build_triangle
Build_triangle(const triSurface &s)
Foam::Vector< scalar >
PolyhedronReader.H
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
Point
CGAL::Point_3< K > Point
Definition: CGALIndexedPolyhedron.H:51
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
Foam::PolyhedronReader::Build_triangle::operator()
void operator()(HDS &hds)