Nektar++
Loading...
Searching...
No Matches
GeomElements.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: GeomElements.cpp
4//
5// For more information, please see: http://www.nektar.info/
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Python wrapper for GeomElements.
32//
33////////////////////////////////////////////////////////////////////////////////
34
37
46
47using namespace Nektar;
48using namespace Nektar::SpatialDomains;
49
51
52template <class T, class S>
53unique_ptr_objpool<T> Geometry_Init(int id, py::list &facets)
54{
55 std::array<S *, T::kNfacets> geomArr;
56
57 if (py::len(facets) != T::kNfacets)
58 {
59 throw new NekError("Incorrect number of facets for this geometry");
60 }
61
62 for (int i = 0; i < T::kNfacets; ++i)
63 {
64 geomArr[i] = py::cast<S *>(facets[i]);
65 }
66
67 return ObjPoolManager<T>::AllocateUniquePtr(id, geomArr);
68}
69
70template <class T, class S>
72 CurveSharedPtr curve)
73{
74 std::array<S *, T::kNfacets> geomArr;
75
76 if (py::len(facets) != T::kNfacets)
77 {
78 throw new NekError("Incorrect number of facets for this geometry");
79 }
80
81 for (int i = 0; i < T::kNfacets; ++i)
82 {
83 geomArr[i] = py::cast<S *>(facets[i]);
84 }
85
86 return ObjPoolManager<T>::AllocateUniquePtr(id, geomArr, curve);
87}
88
89template <class T, class S> void export_Geom_2d(py::module &m, const char *name)
90{
91 py::classh<T, Geometry2D>(m, name)
92 .def(py::init<>(&Geometry_Init<T, S>), py::arg("id"),
93 py::arg("segments") = py::list())
94 .def(py::init<>(&Geometry_Init_Curved<T, S>), py::arg("id"),
95 py::arg("segments"), py::arg("curve"));
96}
97
98template <class T, class S> void export_Geom_3d(py::module &m, const char *name)
99{
100 py::classh<T, Geometry3D>(m, name).def(py::init<>(&Geometry_Init<T, S>),
101 py::arg("id"),
102 py::arg("segments") = py::list());
103}
104
105unique_ptr_objpool<SegGeom> SegGeom_Init(int id, int coordim, py::list &points,
106 CurveSharedPtr curve)
107{
108 std::array<PointGeom *, 2> geomArr;
109
110 if (py::len(points) != 2)
111 {
112 throw new NekError("Incorrect number of facets for this geometry");
113 }
114
115 for (int i = 0; i < 2; ++i)
116 {
117 geomArr[i] = py::cast<PointGeom *>(points[i]);
118 }
119
120 if (!curve)
121 {
122 return ObjPoolManager<SegGeom>::AllocateUniquePtr(id, coordim, geomArr);
123 }
124 else
125 {
126 return ObjPoolManager<SegGeom>::AllocateUniquePtr(id, coordim, geomArr,
127 curve);
128 }
129}
130
132 NekDouble y, NekDouble z)
133{
134 return ObjPoolManager<PointGeom>::AllocateUniquePtr(coordim, id, x, y, z);
135}
136
138{
139 return py::make_tuple(self.x(), self.y(), self.z());
140}
141
142void export_GeomElements(py::module &m)
143{
144 // Geometry dimensioned base classes
145 py::classh<Geometry1D, Geometry>(m, "Geometry1D");
146 py::classh<Geometry2D, Geometry>(m, "Geometry2D")
147 .def("GetCurve", &Geometry2D::GetCurve,
148 py::return_value_policy::reference);
149 py::classh<Geometry3D, Geometry>(m, "Geometry3D");
150
151 // Point geometries
152 py::classh<PointGeom, Geometry>(m, "PointGeom")
153 .def(py::init<>(&PointGeom_Init), py::arg("coordim"), py::arg("id"),
154 py::arg("x"), py::arg("y"), py::arg("z"))
155 .def("GetCoordinates", &PointGeom_GetCoordinates);
156
157 // Segment geometries
158 py::classh<SegGeom, Geometry>(m, "SegGeom")
159 .def(py::init<>(&SegGeom_Init), py::arg("id"), py::arg("coordim"),
160 py::arg("points") = py::list(),
161 py::arg("curve") = CurveSharedPtr())
162 .def("GetCurve", &SegGeom::GetCurve);
163
164 export_Geom_2d<TriGeom, SegGeom>(m, "TriGeom");
165 export_Geom_2d<QuadGeom, SegGeom>(m, "QuadGeom");
166 export_Geom_3d<TetGeom, TriGeom>(m, "TetGeom");
167 export_Geom_3d<PrismGeom, Geometry2D>(m, "PrismGeom");
168 export_Geom_3d<PyrGeom, Geometry2D>(m, "PyrGeom");
169 export_Geom_3d<HexGeom, QuadGeom>(m, "HexGeom");
170}
unique_ptr_objpool< T > Geometry_Init_Curved(int id, py::list &facets, CurveSharedPtr curve)
py::tuple PointGeom_GetCoordinates(const PointGeom &self)
unique_ptr_objpool< SegGeom > SegGeom_Init(int id, int coordim, py::list &points, CurveSharedPtr curve)
void export_Geom_3d(py::module &m, const char *name)
void export_Geom_2d(py::module &m, const char *name)
Nektar::ErrorUtil::NekError NekError
void export_GeomElements(py::module &m)
unique_ptr_objpool< T > Geometry_Init(int id, py::list &facets)
unique_ptr_objpool< PointGeom > PointGeom_Init(int coordim, int id, NekDouble x, NekDouble y, NekDouble z)
boost::call_traits< DataType >::const_reference x() const
Definition NekPoint.hpp:160
boost::call_traits< DataType >::const_reference z() const
Definition NekPoint.hpp:172
boost::call_traits< DataType >::const_reference y() const
Definition NekPoint.hpp:166
static std::unique_ptr< DataType, UniquePtrDeleter > AllocateUniquePtr(const Args &...args)
CurveSharedPtr GetCurve()
Definition SegGeom.h:80
std::shared_ptr< Curve > CurveSharedPtr
Definition Curve.hpp:60
std::unique_ptr< T, typename ObjPoolManager< T >::UniquePtrDeleter > unique_ptr_objpool