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 146 of file Python/StdExpansion.cpp.

147 {
148  py::class_<StdExpansion, std::shared_ptr<StdExpansion>, boost::noncopyable>(
149  "StdExpansion", py::no_init)
150 
151  .def("GetNcoeffs", &StdExpansion::GetNcoeffs)
152  .def("GetTotPoints", &StdExpansion::GetTotPoints)
153  .def("GetBasisType", &StdExpansion::GetBasisType)
154  .def("GetPointsType", &StdExpansion::GetPointsType)
155  .def("GetNverts", &StdExpansion::GetNverts)
156  .def("GetNtraces", &StdExpansion::GetNtraces)
157  .def("DetShapeType", &StdExpansion::DetShapeType)
158  .def("GetShapeDimension", &StdExpansion::GetShapeDimension)
159  .def("Integral", &StdExpansion::Integral)
160 
161  .def("GetBasis", &StdExpansion::GetBasis,
162  py::return_value_policy<py::copy_const_reference>())
163 
164  .def("GenMatrix", &StdExpansion::GenMatrix)
165  .def("GetStdMatrix", &StdExpansion::GetStdMatrix)
166 
167  .def("FwdTrans", &StdExpansion_FwdTrans)
168  .def("BwdTrans", &StdExpansion_BwdTrans)
169  .def("IProductWRTBase", &StdExpansion_IProductWRTBase)
170 
171  .def("PhysEvaluate", &StdExpansion_PhysEvaluate)
172  .def("L2", &StdExpansion_L2)
173  .def("L2", &StdExpansion_L2_Error)
174 
175  .def("GetCoords", &StdExpansion_GetCoords)
176 
177  .def("PhysDeriv", &StdExpansion_PhysDeriv);
178 }
NekDouble StdExpansion_L2_Error(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
py::tuple StdExpansion_PhysDeriv(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &inarray)
Array< OneD, NekDouble > StdExpansion_IProductWRTBase(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
NekDouble StdExpansion_L2(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
py::tuple StdExpansion_GetCoords(StdExpansionSharedPtr exp)
Array< OneD, NekDouble > StdExpansion_BwdTrans(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_FwdTrans(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)

References Nektar::StdRegions::StdExpansion::DetShapeType(), Nektar::StdRegions::StdExpansion::GenMatrix(), Nektar::StdRegions::StdExpansion::GetBasis(), Nektar::StdRegions::StdExpansion::GetBasisType(), Nektar::StdRegions::StdExpansion::GetNcoeffs(), Nektar::StdRegions::StdExpansion::GetNtraces(), 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().

◆ StdExpansion_BwdTrans()

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

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

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

Referenced by export_StdExpansion().

◆ StdExpansion_FwdTrans()

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

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

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

Referenced by export_StdExpansion().

◆ StdExpansion_GetCoords()

py::tuple StdExpansion_GetCoords ( StdExpansionSharedPtr  exp)

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

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

Referenced by export_StdExpansion().

◆ StdExpansion_IProductWRTBase()

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

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

59 {
60  Array<OneD, NekDouble> out(exp->GetNcoeffs());
61  exp->IProductWRTBase(in, out);
62  return out;
63 }

Referenced by export_StdExpansion().

◆ StdExpansion_L2()

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

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

74 {
75  return exp->L2(in);
76 }

Referenced by export_StdExpansion().

◆ StdExpansion_L2_Error()

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

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

81 {
82  return exp->L2(in, err);
83 }

Referenced by export_StdExpansion().

◆ StdExpansion_PhysDeriv()

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

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

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

Referenced by export_StdExpansion().

◆ StdExpansion_PhysEvaluate()

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

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

68 {
69  return exp->PhysEvaluate(coords, physvals);
70 }

Referenced by export_StdExpansion().