Nektar++
Functions
Python/StdExpansion.cpp File Reference
#include <LibUtilities/Python/NekPyConfig.hpp>
#include <StdRegions/StdExpansion.h>

Go to the source code of this file.

Functions

Array< OneD, NekDoubleStdExpansion_FwdTrans (StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
 
Array< OneD, NekDoubleStdExpansion_BwdTrans (StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
 
Array< OneD, NekDoubleStdExpansion_IProductWRTBase (StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
 
NekDouble StdExpansion_PhysEvaluate (StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
 
NekDouble StdExpansion_L2 (StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
 
NekDouble StdExpansion_L2_Error (StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
 
py::tuple StdExpansion_GetCoords (StdExpansionSharedPtr exp)
 
py::tuple StdExpansion_PhysDeriv (StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &inarray)
 
void export_StdExpansion ()
 

Function Documentation

◆ export_StdExpansion()

void export_StdExpansion ( )

Definition at line 152 of file Python/StdExpansion.cpp.

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::StdExpansion::GenMatrix(), Nektar::StdRegions::StdExpansion::GetBasis(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetNedges(), Nektar::StdRegions::StdExpansion::GetNfaces(), Nektar::StdRegions::StdExpansion::GetNverts(), Nektar::StdRegions::StdExpansion::GetPointsType(), Nektar::StdRegions::StdExpansion::GetShapeDimension(), Nektar::StdRegions::StdExpansion::GetStdMatrix(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::StdRegions::StdExpansion::Integral(), StdExpansion_BwdTrans(), StdExpansion_FwdTrans(), StdExpansion_GetCoords(), StdExpansion_IProductWRTBase(), StdExpansion_L2(), StdExpansion_L2_Error(), StdExpansion_PhysDeriv(), and StdExpansion_PhysEvaluate().

Referenced by BOOST_PYTHON_MODULE().

153 {
154  py::class_<StdExpansion,
155  std::shared_ptr<StdExpansion>,
156  boost::noncopyable>(
157  "StdExpansion", py::no_init)
158 
159  .def("GetNcoeffs", &StdExpansion::GetNcoeffs)
160  .def("GetTotPoints", &StdExpansion::GetTotPoints)
161  .def("GetBasisType", &StdExpansion::GetBasisType)
162  .def("GetPointsType", &StdExpansion::GetPointsType)
163  .def("GetNverts", &StdExpansion::GetNverts)
164  .def("GetNedges", &StdExpansion::GetNedges)
165  .def("GetNfaces", &StdExpansion::GetNfaces)
166  .def("DetShapeType", &StdExpansion::DetShapeType)
167  .def("GetShapeDimension", &StdExpansion::GetShapeDimension)
168  .def("Integral", &StdExpansion::Integral)
169 
170  .def("GetBasis", &StdExpansion::GetBasis,
171  py::return_value_policy<py::copy_const_reference>())
172 
173  .def("GenMatrix", &StdExpansion::GenMatrix)
174  .def("GetStdMatrix", &StdExpansion::GetStdMatrix)
175 
176  .def("FwdTrans", &StdExpansion_FwdTrans)
177  .def("BwdTrans", &StdExpansion_BwdTrans)
178  .def("IProductWRTBase", &StdExpansion_IProductWRTBase)
179 
180  .def("PhysEvaluate", &StdExpansion_PhysEvaluate)
181  .def("L2", &StdExpansion_L2)
182  .def("L2", &StdExpansion_L2_Error)
183 
184  .def("GetCoords", &StdExpansion_GetCoords)
185 
186  .def("PhysDeriv", &StdExpansion_PhysDeriv)
187  ;
188 }
Array< OneD, NekDouble > StdExpansion_FwdTrans(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
The base class for all shapes.
Definition: StdExpansion.h:68
NekDouble StdExpansion_L2_Error(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
py::tuple StdExpansion_GetCoords(StdExpansionSharedPtr exp)
py::tuple StdExpansion_PhysDeriv(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &inarray)
Array< OneD, NekDouble > StdExpansion_BwdTrans(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
NekDouble StdExpansion_L2(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
NekDouble StdExpansion_PhysEvaluate(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
Array< OneD, NekDouble > StdExpansion_IProductWRTBase(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)

◆ StdExpansion_BwdTrans()

Array<OneD, NekDouble> StdExpansion_BwdTrans ( StdExpansionSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 50 of file Python/StdExpansion.cpp.

Referenced by export_StdExpansion().

53 {
54  Array<OneD, NekDouble> out(exp->GetTotPoints());
55  exp->BwdTrans(in, out);
56  return out;
57 }

◆ StdExpansion_FwdTrans()

Array<OneD, NekDouble> StdExpansion_FwdTrans ( StdExpansionSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 41 of file Python/StdExpansion.cpp.

Referenced by export_StdExpansion().

44 {
45  Array<OneD, NekDouble> out(exp->GetNcoeffs());
46  exp->FwdTrans(in, out);
47  return out;
48 }

◆ StdExpansion_GetCoords()

py::tuple StdExpansion_GetCoords ( StdExpansionSharedPtr  exp)

Definition at line 91 of file Python/StdExpansion.cpp.

Referenced by export_StdExpansion().

92 {
93  int nPhys = exp->GetTotPoints();
94  int coordim = exp->GetCoordim();
95 
96  std::vector<Array<OneD, NekDouble> > coords(coordim);
97  for (int i = 0; i < coordim; ++i)
98  {
99  coords[i] = Array<OneD, NekDouble>(nPhys);
100  }
101 
102  switch (coordim)
103  {
104  case 1:
105  exp->GetCoords(coords[0]);
106  return py::make_tuple(coords[0]);
107  break;
108  case 2:
109  exp->GetCoords(coords[0], coords[1]);
110  return py::make_tuple(coords[0], coords[1]);
111  break;
112  case 3:
113  exp->GetCoords(coords[0], coords[1], coords[2]);
114  return py::make_tuple(coords[0], coords[1], coords[2]);
115  break;
116  }
117 
118  return py::tuple();
119 }

◆ StdExpansion_IProductWRTBase()

Array<OneD, NekDouble> StdExpansion_IProductWRTBase ( StdExpansionSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 59 of file Python/StdExpansion.cpp.

Referenced by export_StdExpansion().

62 {
63  Array<OneD, NekDouble> out(exp->GetNcoeffs());
64  exp->BwdTrans(in, out);
65  return out;
66 }

◆ StdExpansion_L2()

NekDouble StdExpansion_L2 ( StdExpansionSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 76 of file Python/StdExpansion.cpp.

Referenced by export_StdExpansion().

79 {
80  return exp->L2(in);
81 }

◆ StdExpansion_L2_Error()

NekDouble StdExpansion_L2_Error ( StdExpansionSharedPtr  exp,
const Array< OneD, const NekDouble > &  in,
const Array< OneD, const NekDouble > &  err 
)

Definition at line 83 of file Python/StdExpansion.cpp.

Referenced by export_StdExpansion().

87 {
88  return exp->L2(in, err);
89 }

◆ StdExpansion_PhysDeriv()

py::tuple StdExpansion_PhysDeriv ( StdExpansionSharedPtr  exp,
const Array< OneD, const NekDouble > &  inarray 
)

Definition at line 121 of file Python/StdExpansion.cpp.

Referenced by export_StdExpansion().

123 {
124  int nPhys = exp->GetTotPoints();
125  int coordim = exp->GetCoordim();
126 
127  std::vector<Array<OneD, NekDouble> > derivs(coordim);
128  for (int i = 0; i < coordim; ++i)
129  {
130  derivs[i] = Array<OneD, NekDouble>(nPhys);
131  }
132 
133  switch (coordim)
134  {
135  case 1:
136  exp->PhysDeriv(inarray, derivs[0]);
137  return py::make_tuple(derivs[0]);
138  break;
139  case 2:
140  exp->PhysDeriv(inarray, derivs[0], derivs[1]);
141  return py::make_tuple(derivs[0], derivs[1]);
142  break;
143  case 3:
144  exp->PhysDeriv(inarray, derivs[0], derivs[1], derivs[2]);
145  return py::make_tuple(derivs[0], derivs[1], derivs[2]);
146  break;
147  }
148 
149  return py::tuple();
150 }

◆ StdExpansion_PhysEvaluate()

NekDouble StdExpansion_PhysEvaluate ( StdExpansionSharedPtr  exp,
const Array< OneD, const NekDouble > &  coords,
const Array< OneD, const NekDouble > &  physvals 
)

Definition at line 68 of file Python/StdExpansion.cpp.

Referenced by export_StdExpansion().

72 {
73  return exp->PhysEvaluate(coords, physvals);
74 }