COMPASS  5.0.0
End-to-end AO simulation tool using GPU acceleration
CarmaObjInterfacer Struct Reference

#include <obj.hpp>

Collaboration diagram for CarmaObjInterfacer:

Static Public Member Functions

template<typename T >
static void call (py::module &mod)
 

Detailed Description

Definition at line 59 of file obj.hpp.

Member Function Documentation

◆ call()

template<typename T >
static void CarmaObjInterfacer::call ( py::module &  mod)
static

< sum

Definition at line 61 of file obj.hpp.

61  {
62  auto name = appendName<T>("obj_");
63  using Class = CarmaObj<T>;
64  using ClassHost = CarmaHostObj<T>;
65 
66  py::class_<Class>(mod, name.data(), py::buffer_protocol())
67  .def(py::init([](CarmaContext &c,
68  const py::array_t<T, py::array::f_style |
69  py::array::forcecast> &data) {
70  int ndim = data.ndim() + 1;
71  std::vector<long> data_dims(ndim);
72  data_dims[0] = data.ndim();
73  copy(data.shape(), data.shape() + data.ndim(),
74  begin(data_dims) + 1);
75  return std::unique_ptr<Class>(
76  new Class(&c, data_dims.data(), (const T *)data.data()));
77  }),
78  "TODO", // TODO do the documentation...
79  py::arg("context").none(false), py::arg("h_data").none(false))
80 
81  .def(py::init([](CarmaContext &c, const Class &data) {
82  return std::unique_ptr<Class>(new Class(&c, &data));
83  }),
84  "TODO", // TODO do the documentation...
85  py::arg("context").none(false), py::arg("d_data").none(false))
86 
87  .def_buffer([](Class &frame) -> py::buffer_info {
88  frame.sync_h_data();
89 
90  const long *dims = frame.get_dims();
91  std::vector<ssize_t> shape(dims[0]);
92  std::vector<ssize_t> strides(dims[0]);
93  ssize_t stride = sizeof(T);
94 
95  // C-style
96  // for (ssize_t dim(dims[0] - 1); dim >= 0; --dim)
97  // {
98  // shape[dim] = dims[dim + 1];
99  // strides[dim] = stride;
100  // stride *= shape[dim];
101  // }
102 
103  // F-style
104  for (ssize_t dim(0); dim < dims[0]; ++dim) {
105  shape[dim] = dims[dim + 1];
106  strides[dim] = stride;
107  stride *= shape[dim];
108  }
109 
110  return py::buffer_info(frame.get_h_data(), sizeof(T),
111  py::format_descriptor<T>::format(), dims[0],
112  shape, strides);
113  })
114 
115  .def("__repr__", &Class::to_string)
116 
117  // int get_nb_streams()
118  .def_property_readonly("nb_streams", &Class::get_nb_streams,
119  "TODO") // TODO do the documentation...
120  // int add_stream()
121  .def("add_stream", (int (Class::*)()) & Class::add_stream,
122  "TODO") // TODO do the documentation...
123  // int add_stream(int nb)
124  .def("add_stream", (int (Class::*)(int)) & Class::add_stream, "TODO",
125  py::arg("np")) // TODO do the documentation...
126  // int del_stream()
127  .def("del_stream", (int (Class::*)()) & Class::del_stream,
128  "TODO") // TODO do the documentation...
129  // int del_stream(int nb)
130  .def("del_stream", (int (Class::*)(int)) & Class::del_stream, "TODO",
131  py::arg("np")) // TODO do the documentation...
132  // int wait_stream(int stream)
133  .def("wait_stream", &Class::wait_stream, "TODO",
134  py::arg("steam")) // TODO do the documentation...
135  .def("swap_ptr", [](Class &obj, Class &obj2){
136  obj.swap_ptr(obj2.get_data());
137  }, "TODO",
138  py::arg("ptr")) // TODO do the documentation...
139  // int wait_all_streams()
140 #ifdef USE_OCTOPUS
141  .def("swap_ptr", [](Class &obj, ipc::Cacao<T> &obj2){
142  obj.swap_ptr(obj2.inputPtr());
143  }, "TODO",
144  py::arg("ptr")) // TODO do the documentation...
145  // int wait_all_streams()
146 #endif
147  .def("wait_all_streams", &Class::wait_all_streams,
148  "TODO") // TODO do the documentation...
149 
150  // const long *get_dims()
151  .def_property_readonly("shape",
152  [](Class &frame) -> py::array_t<long> {
153  long nb_dim = frame.get_dims(0);
154  const long *c_dim = frame.get_dims() + 1;
155  return py::array_t<long>(nb_dim, c_dim);
156  },
157  "TODO") // TODO do the documentation...
158 
159  // int get_nb_elements()
160  .def_property_readonly("nbElem", &Class::get_nb_elements,
161  "TODO") // TODO do the documentation...
162  // CarmaContext* get_context()
163  .def_property_readonly("context", &Class::get_context,
164  "TODO") // TODO do the documentation...
165  // int get_device()
166  .def_property_readonly("device", &Class::get_device,
167  "TODO") // TODO do the documentation...
168  // int get_o_data()
169  .def_property_readonly("o_data", &Class::get_o_data_value,
170  "TODO") // TODO do the documentation...
171 
172  // int host2device(T_data *data);
173  .def("host2device",
174  [](Class &c,
175  py::array_t<T, py::array::f_style | py::array::forcecast>
176  &data) { c.host2device((const T *)data.data()); },
177  "TODO",
178  py::arg("data").none(false)) // TODO do the documentation...
179  // int device2host(T_data *data);
180  .def("device2host",
181  [](Class &c,
182  py::array_t<T, py::array::f_style | py::array::forcecast>
183  &data) { c.device2host((T *)data.mutable_data()); },
184  "TODO",
185  py::arg("data").none(false)) // TODO do the documentation...
186 
187  // int copy_into(T_data *data, int nb_elem);
188  .def("copy_into",
189  [](Class &src, Class &dest, long nb_elem) {
190  if (nb_elem < 0) {
191  nb_elem = src.get_nb_elements();
192  }
193  src.copy_into(dest, nb_elem);
194  },
195  "TODO", py::arg("dest"),
196  py::arg("nb_elem") = -1) // TODO do the documentation...
197  // int copy_from(T_data *data, int nb_elem);
198  .def("copy_from",
199  [](Class &dest, Class &src, long nb_elem) {
200  if (nb_elem < 0) {
201  nb_elem = dest.get_nb_elements();
202  }
203  dest.copy_from(src, nb_elem);
204  },
205  "TODO", py::arg("data"),
206  py::arg("nb_elem") = -1) // TODO do the documentation...
207 #ifdef USE_OCTOPUS
208  .def("copy_into",(int (Class::*)(ipc::Cacao<T>*))&Class::copy_into)
209  .def("copy_from",(int (Class::*)(ipc::Cacao<T>*))&Class::copy_from)
210 #endif
211  // inline int reset()
212  .def("reset", &Class::reset, "TODO") // TODO do the documentation...
213 
215  // T_data sum();
216  .def("sum", &Class::sum, "TODO") // TODO do the documentation...
217  // void init_reduceCub();
218  .def("init_reduceCub", &Class::init_reduceCub,
219  "TODO") // TODO do the documentation...
220  // void reduceCub();
221  .def("reduceCub", &Class::reduceCub,
222  "TODO") // TODO do the documentation...
223  // void clip(T_data min, T_data max);
224  .def("clip", &Class::clip, "TODO", py::arg("data_min").none(false),
225  py::arg("data_max").none(false)) // TODO do the documentation...
226 
227  // /**< transpose */
228  // int transpose(CarmaObj<T_data> *source);
229  .def("transpose", &Class::transpose, "TODO",
230  py::arg("source").none(false)) // TODO do the documentation...
231  // //CarmaObj<T_data>& operator= (const CarmaObj<T_data>& obj);
232 
233  // /**< Cublas V2 */
234  // int imax(int incx);
235  .def("aimax", &Class::aimax, "TODO",
236  py::arg("incx") = 1) // TODO do the documentation...
237  // int imin(int incx);
238  .def("aimin", &Class::aimin, "TODO",
239  py::arg("incx") = 1) // TODO do the documentation...
240  // T_data asum(int incx);
241  .def("asum", &Class::asum, "TODO",
242  py::arg("incx") = 1) // TODO do the documentation...
243  // T_data nrm2(int incx);
244  .def("nrm2", &Class::nrm2, "TODO",
245  py::arg("incx") = 1) // TODO do the documentation...
246  // T_data dot(CarmaObj<T_data> *source, int incx, int incy);
247  .def("dot", &Class::dot, "TODO", py::arg("source").none(false),
248  py::arg("incx") = 1,
249  py::arg("incy") = 1) // TODO do the documentation...
250  // void scale(T_data alpha, int incx);
251  .def("scale", &Class::scale, "TODO", py::arg("scale").none(false),
252  py::arg("incx") = 1) // TODO do the documentation...
253  // void swap(CarmaObj<T_data> *source, int incx, int incy);
254  .def("swap", &Class::swap, "TODO", py::arg("source").none(false),
255  py::arg("incx") = 1,
256  py::arg("incy") = 1) // TODO do the documentation...
257  // void copy(CarmaObj<T_data> *source, int incx, int incy);
258  .def("copy", &Class::copy, "TODO") // TODO do the documentation...
259  // void axpy(T_data alpha, CarmaObj<T_data> *source, int incx, int
260  // incy);
261  .def("axpy", &Class::axpy, "TODO", py::arg("alpha"),
262  py::arg("source").none(false), py::arg("incx") = 1,
263  py::arg("incy") = 1, py::arg("offset") = 0) // TODO do the documentation...
264  // void rot(CarmaObj<T_data> *source, int incx, int incy, T_data sc,
265  // T_data ss);
266  .def("rot", &Class::rot, "TODO") // TODO do the documentation...
267 
268  // void gemv(char trans, T_data alpha, CarmaObj<T_data> *matA, int lda,
269  // CarmaObj<T_data> *vectx, int incx, T_data beta, int incy);
270  .def("magma_gemv",
271  [](Class &mat, Class &vectx, T alpha, char op, Class *vecty,
272  T beta) {
273  if (vecty == nullptr) {
274  long dims[] = {1, 0};
275  if (op == 'N' || op == 'n') {
276  dims[1] = mat.get_dims(1);
277  } else {
278  dims[1] = mat.get_dims(2);
279  }
280  vecty = new Class(mat.get_context(), dims);
281  vecty->reset();
282  }
283  carma_magma_gemv(op, mat.get_dims(1), mat.get_dims(2), alpha, mat.get_data(), mat.get_dims(1), vectx.get_data(), 1, beta, vecty->get_data(), 1);
284  return vecty;
285  },
286  "this method performs one of the matrix‐vector operations vecty = "
287  "alpha * op(mat) * vectx + beta * vecty",
288  py::arg("vectx"), py::arg("alpha") = 1, py::arg("op") = 'N',
289  py::arg("vecty") = nullptr, py::arg("beta") = 0) // &Class::gemv)
290  .def("gemv",
291  [](Class &mat, Class &vectx, T alpha, char op, Class *vecty,
292  T beta) {
293  if (vecty == nullptr) {
294  long dims[] = {1, 0};
295  if (op == 'N' || op == 'n') {
296  dims[1] = mat.get_dims(1);
297  } else {
298  dims[1] = mat.get_dims(2);
299  }
300  vecty = new Class(mat.get_context(), dims);
301  vecty->reset();
302  }
303  vecty->gemv(op, alpha, &mat, mat.get_dims(1), &vectx, 1, beta, 1);
304  return vecty;
305  },
306  "this method performs one of the matrix‐vector operations vecty = "
307  "alpha * op(mat) * vectx + beta * vecty",
308  py::arg("vectx"), py::arg("alpha") = 1, py::arg("op") = 'N',
309  py::arg("vecty") = nullptr, py::arg("beta") = 0) // &Class::gemv)
310  // void ger(T_data alpha, CarmaObj<T_data> *vectx, int incx,
311  // CarmaObj<T_data> *vecty, int incy, int lda);
312 #ifdef CAN_DO_HALF
313  .def("gemv",
314  [](CarmaObj<half> &mat, CarmaObj<half> &vectx, float alpha,
315  char op, CarmaObj<half> *vecty, float beta) {
316  if (vecty == nullptr) {
317  long dims[] = {1, 0, 1};
318  if (op == 'N' || op == 'n') {
319  dims[1] = mat.get_dims(1);
320  } else {
321  dims[1] = mat.get_dims(2);
322  }
323  vecty = new CarmaObj<half>(mat.get_context(), dims);
324  vecty->reset();
325  }
326  vecty->gemv(op, __float2half(alpha), &mat, mat.get_dims(1),
327  &vectx, 1, __float2half(beta), 1);
328  return vecty;
329  },
330  "this method performs one of the matrix‐vector operations vecty = "
331  "alpha * op(mat) * vectx + beta * vecty",
332  py::arg("vectx"), py::arg("alpha") = 1, py::arg("op") = 'N',
333  py::arg("vecty") = nullptr, py::arg("beta") = 0) // &Class::gemv)
334 #endif
335  // void ger(T_data alpha, CarmaObj<T_data> *vectx, int incx,
336  // CarmaObj<T_data> *vecty, int incy, int lda);
337  .def("ger",
338  [](Class &vectx, Class &vecty, Class *mat, T alpha) {
339  std::unique_ptr<Class> ptr_res;
340  if (mat == nullptr) {
341  long dims[] = {2, vectx.get_nb_elements(), vecty.get_nb_elements()};
342  mat = new Class(vectx.get_context(), dims);
343  mat->reset();
344  }
345  mat->ger(alpha, &vectx, 1, &vecty, 1, vectx.get_nb_elements());
346  return mat;
347  },
348  "this method performs the symmetric rank 1 operation A = alpha * "
349  "x * y T + A",
350  py::arg("vecty"), py::arg("mat") = nullptr,
351  py::arg("alpha") = 1) // &Class::ger)
352  // void symv(char uplo, T_data alpha, CarmaObj<T_data> *matA,
353  // int lda, CarmaObj<T_data> *vectx, int incx, T_data beta,
354  // int incy);
355  .def("symv",
356  [](Class &mat, Class &vectx, T alpha, char uplo, Class *vecty,
357  T beta) {
358  int lda = mat.get_dims(2);
359  if (vecty == nullptr) {
360  long dims[] = {1, lda};
361  vecty = new Class(mat.get_context(), dims);
362  vecty->reset();
363  }
364  vecty->symv(uplo, alpha, &mat, lda, &vectx, 1, beta, 1);
365  return vecty;
366  },
367  "this method performs one of the matrix‐vector operations vecty = "
368  "alpha * mat * vectx + beta * vecty",
369  py::arg("vectx"), py::arg("alpha") = 1, py::arg("uplo") = 'l',
370  py::arg("vecty") = nullptr, py::arg("beta") = 0) // &Class::gemv)
371  // void gemm(char transa, char transb, T_data alpha, CarmaObj<T_data>
372  // *matA,
373  // int lda, CarmaObj<T_data> *matB, int ldb, T_data beta, int
374  // ldc);
375  .def("gemm",
376  [](Class &matA, Class &matB, char op_a, char op_b, T alpha,
377  Class *matC, T beta) /*-> std::unique_ptr<Class>*/ {
378  int lda, ldb, ldc, m, n, k;
379  if (op_a == 'N' || op_a == 'n') {
380  m = matA.get_dims(1);
381  k = matA.get_dims(2);
382  } else {
383  m = matA.get_dims(2);
384  k = matA.get_dims(1);
385  }
386 
387  if (op_b == 'N' || op_b == 'n') {
388  k = matB.get_dims(1);
389  n = matB.get_dims(2);
390  } else {
391  k = matB.get_dims(2);
392  n = matB.get_dims(1);
393  }
394 
395  if (matC == nullptr) {
396  long dims[] = {2, m, n};
397  matC = new Class(matA.get_context(), dims);
398  matC->reset();
399  }
400 
401  carma_gemm<T>(matA.get_context()->get_cublas_handle(), op_a, op_b,
402  m, n, k, alpha, matA, matA.get_dims(1), matB,
403  matB.get_dims(1), beta, *matC, matC->get_dims(1));
404  return matC;
405  },
406  "this method performs one of the matrix‐marix operations matC = "
407  "alpha * op_a(matA) * op_b(matB) + beta * matC",
408  py::arg("matB"), py::arg("op_a") = 'N', py::arg("op_b") = "N",
409  py::arg("alpha") = 1, py::arg("matC") = nullptr,
410  py::arg("beta") = 0)
411 
412  // void symm(char side, char uplo, T_data alpha,
413  // CarmaObj<T_data> *matA, int lda, CarmaObj<T_data> *matB,
414  // int ldb, T_data beta, int ldc);
415  .def("symm",
416  [](Class &matA, Class &matB, T alpha, Class *matC, T beta,
417  char side, char uplo) {
418  if (matC == nullptr) {
419  long dims[] = {2, matB.get_dims(1), matB.get_dims(2)};
420  matC = new Class(matA.get_context(), dims);
421  }
422  carma_symm<T>(matA.get_context()->get_cublas_handle(), side, uplo,
423  matB.get_dims(1), matB.get_dims(2), alpha, matA,
424  matA.get_dims(1), matB, matB.get_dims(1), beta,
425  *matC, matC->get_dims(1));
426  return matC;
427  // matA.symm(uplo, alpha, &matB, lda, &vectx, 1, beta, 1);
428  },
429  "this method performs one of the matrix‐marix operations matC = "
430  "alpha * matA * matB + beta * C",
431  py::arg("matB"), py::arg("alpha") = 1, py::arg("matC") = nullptr,
432  py::arg("beta") = 0, py::arg("side") = "l", py::arg("uplo") = "u")
433 
434  /*
435  template <class T_data>
436  cublasStatus_t carma_symm(cublasHandle_t cublas_handle, char side,
437  char uplo, int m, int n, T_data alpha,
438  T_data *matA, int lda, T_data *matB, int ldb,
439  T_data beta, T_data *matC, int ldc);
440  */
441 
442  // void syrk(char uplo, char transa, T_data alpha,
443  // CarmaObj<T_data> *matA, int lda, T_data beta, int ldc);
444  .def("syrk",
445  [](Class &matA, char fill, char op, T alpha, Class *matC, T beta) {
446  int n, k;
447  if (op == 'N' || op == 'n') {
448  n = matA.get_dims(1);
449  k = matA.get_dims(2);
450  } else {
451  n = matA.get_dims(2);
452  k = matA.get_dims(1);
453  }
454  if (matC == nullptr) {
455  long dims[] = {2, n, n};
456  matC = new Class(matA.get_context(), dims);
457  matC->reset();
458  }
459  carma_syrk<T>(matA.get_context()->get_cublas_handle(), fill, op, n,
460  k, alpha, matA, matA.get_dims(1), beta, *matC,
461  matC->get_dims(1));
462  return matC;
463  },
464  "this method performs the symmetric rank- k update",
465  py::arg("fill") = "U", py::arg("op") = 'N', py::arg("alpha") = 1,
466  py::arg("matC") = nullptr, py::arg("beta") = 0)
467  // void syrkx(char uplo, char transa, T_data alpha,
468  // CarmaObj<T_data> *matA, int lda, CarmaObj<T_data> *matB,
469  // int ldb, T_data beta, int ldc);
470  .def("syrkx",
471  [](Class &matA, Class &matB, char fill, char op, T alpha,
472  Class *matC, T beta) {
473  int n, k;
474  if (op == 'N' || op == 'n') {
475  n = matA.get_dims(1);
476  k = matA.get_dims(2);
477  } else {
478  n = matA.get_dims(2);
479  k = matA.get_dims(1);
480  }
481  if (matC == nullptr) {
482  long dims[] = {2, n, n};
483  matC = new Class(matA.get_context(), dims);
484  matC->reset();
485  }
486  carma_syrkx<T>(matA.get_context()->get_cublas_handle(), fill, op,
487  n, k, alpha, matA, matA.get_dims(1), matB,
488  matB.get_dims(1), beta, *matC, matC->get_dims(1));
489  return matC;
490  },
491  "this method performs the symmetric rank- k update",
492  py::arg("matB"), py::arg("fill") = "U", py::arg("op") = 'N',
493  py::arg("alpha") = 1, py::arg("matC") = nullptr,
494  py::arg("beta") = 0)
495  // void geam(char transa, char transb, T_data alpha, CarmaObj<T_data>
496  // *matA,
497  // int lda, T_data beta, CarmaObj<T_data> *matB, int ldb, int
498  // ldc);
499  .def("geam",
500  [](Class &matA, Class &matB, char opA, char opB, T alpha,
501  Class *matC, T beta) {
502  int m, n;
503  if (opA == 'N' || opA == 'n') {
504  m = matA.get_dims(1);
505  n = matA.get_dims(2);
506  } else {
507  m = matA.get_dims(2);
508  n = matA.get_dims(1);
509  }
510  if (matC == nullptr) {
511  long dims[] = {2, m, n};
512  matC = new Class(matA.get_context(), dims);
513  matC->reset();
514  }
515  carma_geam<T>(matA.get_context()->get_cublas_handle(), opA, opB, m,
516  n, alpha, matA, matA.get_dims(1), beta, matB,
517  matB.get_dims(1), *matC, matC->get_dims(1));
518  return matC;
519  },
520  "this method performs the symmetric rank- k update",
521  py::arg("matB"), py::arg("opA") = 'N', py::arg("opB") = 'N',
522  py::arg("alpha") = 1, py::arg("matC") = nullptr,
523  py::arg("beta") = 0)
524  .def("dgmm",
525  [](Class &matA, Class &vectX, T alpha, char side, Class *matC,
526  int incx) {
527  if (matC == nullptr) {
528  long dims[] = {2, matA.get_dims(1), matA.get_dims(2)};
529  matC = new Class(matA.get_context(), dims);
530  matC->reset();
531  }
532  carma_dgmm<T>(matA.get_context()->get_cublas_handle(), side,
533  matA.get_dims(1), matA.get_dims(2), matA,
534  matA.get_dims(1), vectX, incx, *matC,
535  matC->get_dims(1));
536  return matC;
537  },
538  "this method performs one of the matrix‐marix operations matC = "
539  "diag(vectX)*matA if side='l'",
540  py::arg("vectX"), py::arg("alpha") = 1, py::arg("side") = "r",
541  py::arg("matC") = nullptr, py::arg("incx") = 1)
542 
543  // /**< Curand */
544  .def("is_rng_init", &Class::is_rng_init)
545  // int init_prng();
546  .def("init_prng", (int (Class::*)()) & Class::init_prng)
547  // int init_prng(long seed);
548  .def("init_prng", (int (Class::*)(long)) & Class::init_prng)
549  // int destroy_prng();
550  .def("destroy_prng", &Class::destroy_prng)
551  // int prng(T_data *output, char gtype, float alpha, float beta);
552  .def("prng", (int (Class::*)(T *, char, float, float)) & Class::prng)
553  // int prng(T_data *output, char gtype, float alpha);
554  .def("prng", (int (Class::*)(T *, char, float)) & Class::prng)
555  // int prng(char gtype, float alpha, float beta);
556  .def("prng", (int (Class::*)(char, float, float)) & Class::prng)
557  // int prng(char gtype, float alpha);
558  .def("prng", (int (Class::*)(char, float)) & Class::prng)
559  // int prng(char gtype);
560  .def("prng", (int (Class::*)(char)) & Class::prng)
561 
562  .def("random",
563  [](Class &data, int seed, char gtype) {
564  data.init_prng(seed);
565  data.prng(gtype);
566  },
567  py::arg("seed") = 1234, py::arg("j") = 'U')
568 
569  .def("random_host",
570  [](Class &data, int seed, char gtype) {
571  data.init_prng_host(seed);
572  data.prng_host(gtype);
573  },
574  py::arg("seed") = 1234, py::arg("j") = 'U')
575 
576  // int prng_montagn( float init_montagn );
577  .def("prng_montagn", &Class::prng_montagn)
578 
579  // int init_prng_host(int seed);
580  .def("init_prng_host", (int (Class::*)(int)) & Class::init_prng_host)
581  // int prng_host(char gtype);
582  .def("prng_host", (int (Class::*)(char)) & Class::prng_host)
583  // int prng_host(char gtype, T_data stddev);
584  .def("prng_host", (int (Class::*)(char, T)) & Class::prng_host)
585  // int prng_host(char gtype, T_data stddev, T_data alpha);
586  .def("prng_host", (int (Class::*)(char, T, T)) & Class::prng_host)
587  // int destroy_prng_host();
588  .def("destroy_prng_host", &Class::destroy_prng_host)
589 
590  .def("fft",
591  [](Class &data, Class &dest, int direction) {
592  throw std::runtime_error("not implemented");
593  // const long *dims = data.get_dims();
594  // cufftHandle *handle = data.get_plan();
595  // if(dest == nullptr) {
596  // dest = Class(data.get_context(), dims);
597  // }
598  // carma_initfft(dims, handle, carma_select_plan<T,T>());
599  // CarmaFFT(data.get_data(), dest.get_data(), direction, handle);
600  },
601  py::arg("dest") = nullptr, py::arg("direction") = 1);
602  // CU functions clip
603  // template<class T_data>
604  // void clip_array(T_data *d_data, T_data min, T_data max, int N,
605  // CarmaDevice *device);
606 
607  // CU functions sum
608  // template<class T_data>
609  // void reduce(int size, int threads, int blocks, T_data *d_idata,
610  // T_data *d_odata);
611  // template<class T_data>
612  // T_data reduce(T_data * data, int N);
613 
614  // CU functions transpose
615  // template<class T_data>
616  // int transposeCU(T_data *d_idata, T_data *d_odata, long N1, long N2);
617 
618  // CU functions generic
619  // template<class T_data>
620  // int launch_generic1d(T_data *d_idata, T_data *d_odata, int N,
621  // CarmaDevice *device);
622  // template<class T_data>
623  // int launch_generic2d(T_data *d_odata, T_data *d_idata, int N1, int N2);
624 
625  // CU functions curand
626  // int carma_prng_init(int *seed, const int nb_threads, const int nb_blocks,
627  // curandState *state);
628  // template<class T>
629  // int carma_prng_cu(T *results, const int nb_threads, const int nb_blocks,
630  // curandState *state, char gtype, int n, float alpha,
631  // float beta);
632  // template<class T>
633  // int carma_curand_montagn(curandState *state, T *d_odata, int N,
634  // CarmaDevice *device);
635 
636  // CU functions fft
637  // template<class T_in, class T_out>
638  // cufftType carma_select_plan();
639  // template<class T_in, class T_out>
640  // void carma_initfft(const long *dims_data, cufftHandle *plan, cufftType
641  // type_plan); template<class T_in, class T_out> int CarmaFFT(T_in *input,
642  // T_out *output, int dir, cufftHandle plan);
643 
644  // CU functions generic
645  // template<class T_data>
646  // int fillindex(T_data *d_odata, T_data *d_idata, int *indx, int N,
647  // CarmaDevice *device);
648  // template<class T_data>
649  // int fillvalues(T_data *d_odata, T_data val, int N,
650  // CarmaDevice *device);
651  // template<class T>
652  // int getarray2d(T *d_odata, T *d_idata, int x0, int Ncol, int NC, int N,
653  // CarmaDevice *device);
654  // template<class T>
655  // int fillarray2d(T *d_odata, T *d_idata, int x0, int Ncol, int NC, int N,
656  // CarmaDevice *device);
657  // template<class T>
658  // int fillarray2d2(T *d_odata, T *d_idata, int x0, int Ncol, int NC, int N,
659  // CarmaDevice *device);
660  // template<class T>
661  // int fill_sym_matrix(char src_uplo, T *d_data, int Ncol, int N,
662  // CarmaDevice *device);
663  // template<class T>
664  // int carma_plus(T *d_odata, T elpha, int N, CarmaDevice *device);
665  // template<class T>
666  // int carma_plusai(T *d_odata, T *i_data, int i, int sgn, int N,
667  // CarmaDevice *device);
668 
669  // CU functions fftconv
670  // int fftconv_unpad(float *d_odata, float *d_idata, int fftW, int dataH,
671  // int dataW, int N, int n, int nim);
672  // int carma_initfftconv(CarmaObjS *data_in, CarmaObjS *kernel_in, CarmaObjS
673  // *padded_data, CarmaObjC *padded_spectrum, int kernelY, int kernelX);
674 
675  // CPP functions fftconv
676  // int carma_fftconv(CarmaObjS *data_out, CarmaObjS *padded_data,
677  // CarmaObjC *padded_spectrum, int kernelY, int kernelX);
678 
679  // MAGMA functions
680 
681  // template<class T>
682  // int carma_svd(CarmaObj<T> *imat, CarmaObj<T> *eigenvals,
683  // CarmaObj<T> *mod2act, CarmaObj<T> *mes2mod);
684  // mod.def(appendName<T>("svd_").data(), &carma_svd<T>);
685 
686  // TODO after CarmaHostObj
687  // template<class T>
688  // int carma_magma_syevd(char jobz, CarmaObj<T> *mat, CarmaHostObj<T>
689  // *eigenvals);
690  // mod.def( appendName<T>("syevd_").data(), py::overload_cast<char, Class *,
691  // ClassHost *>(&carma_magma_syevd<T>)); mod.def(
692  // appendName<T>("syevd_").data(), py::overload_cast<char, long, T *, T
693  // *>(&carma_magma_syevd<T>));
694  mod.def(
695  appendName<T>("magma_syevd_").data(),
696  [](Class &d_mat_a, ClassHost &eigenvals, Class *d_U, bool computeU) {
697  if (d_U == nullptr) {
698  if (computeU) {
699  carma_magma_syevd('V', &d_mat_a, &eigenvals);
700  } else {
701  carma_magma_syevd('N', &d_mat_a, &eigenvals);
702  }
703  } else {
704  d_U->copy_from(d_mat_a, d_mat_a.get_nb_elements());
705  if (computeU) {
706  carma_magma_syevd('V', d_U, &eigenvals);
707  } else {
708  carma_magma_syevd('N', d_U, &eigenvals);
709  }
710  }
711  },
712  py::arg("d_mat_a"), py::arg("eigenvals"), py::arg("d_U") = nullptr,
713  py::arg("computeU") = true);
714 
715  mod.def(
716  appendName<T>("syevd_").data(),
717  [](Class &d_mat_a, Class &eigenvals, Class *d_U, bool computeU) {
718  if (d_U == nullptr) {
719  if (computeU) {
720  carma_syevd(CUSOLVER_EIG_MODE_VECTOR, &d_mat_a, &eigenvals);
721  } else {
722  carma_syevd(CUSOLVER_EIG_MODE_NOVECTOR, &d_mat_a, &eigenvals);
723  }
724  } else {
725  d_U->copy_from(d_mat_a, d_mat_a.get_nb_elements());
726  if (computeU) {
727  carma_syevd(CUSOLVER_EIG_MODE_VECTOR, d_U, &eigenvals);
728  } else {
729  carma_syevd(CUSOLVER_EIG_MODE_NOVECTOR, d_U, &eigenvals);
730  }
731  }
732  },
733  py::arg("d_mat_a"), py::arg("eigenvals"), py::arg("d_U") = nullptr,
734  py::arg("computeU") = true);
735  // template<class T, int method>
736  // int carma_magma_syevd(char jobz, CarmaObj<T> *mat, CarmaHostObj<T>
737  // *eigenvals); template<class T> int carma_magma_syevd_m(long ngpu, char
738  // jobz, long N, T *mat, T *eigenvals); template<class T> int
739  // carma_magma_syevd_m(long ngpu, char jobz, CarmaHostObj<T> *mat,
740  // CarmaHostObj<T> *eigenvals);
741  // template<class T>
742  // int carma_magma_syevd_m(long ngpu, char jobz, CarmaHostObj<T> *mat,
743  // CarmaHostObj<T> *eigenvals, CarmaHostObj<T> *U);
744  // template<class T>
745  // int carma_magma_getri(CarmaObj<T> *d_iA);
746  mod.def(appendName<T>("magma_getri_").data(), &carma_magma_getri<T>);
747 
748  // template<class T>
749  // int carma_magma_potri(CarmaObj<T> *d_iA);
750  mod.def(appendName<T>("magma_potri_").data(), &carma_magma_potri<T>);
751 
752  // TODO after CarmaHostObj
753  // template<class T>
754  // int carma_magma_potri_m(long num_gpus, CarmaHostObj<T> *h_A,
755  // CarmaObj<T> *d_iA);
756 
757  // MAGMA functions (direct access)
758  // template<class T>
759  // int carma_magma_syevd(char jobz, long N, T *mat, T *eigenvals);
760  // template<class T, int method>
761  // int carma_magma_syevd(char jobz, long N, T *mat, T *eigenvals);
762  // template<class T>
763  // int carma_magma_syevd_m(long ngpu, char jobz, long N, T *mat, T
764  // *eigenvals); template<class T> int carma_magma_potri_m(long num_gpus,
765  // long N, T *h_A, T *d_iA);
766 
767  }
Here is the call graph for this function:

The documentation for this struct was generated from the following file:
carma_magma_gemv
int carma_magma_gemv(char trans, int m, int n, T_data alpha, T_data *matA, int lda, T_data *vectx, int incx, T_data beta, T_data *vecty, int incy)
CarmaObj::get_dims
const long * get_dims()
Definition: carma_obj.h:239
test_cusolver.n
int n
Definition: test_cusolver.py:7
CarmaObj::gemv
void gemv(char trans, T_data alpha, CarmaObj< T_data > *matA, int lda, CarmaObj< T_data > *vectx, int incx, T_data beta, int incy)
CarmaObj::get_context
CarmaContext * get_context()
Definition: carma_obj.h:242
correlation_bokeh.data
data
Definition: correlation_bokeh.py:280
test_cublas1.c
c
Definition: test_cublas1.py:14
carma_syevd
int carma_syevd(cusolverEigMode_t jobz, CarmaObj< T > *mat, CarmaObj< T > *eigenvals)
carma_utils::to_string
std::string to_string(const T &n)
Definition: carma_utils.h:82
test_cusolver.m
int m
Definition: test_cusolver.py:6
CarmaObj< T >
CarmaContext
this class provides the context in which CarmaObj are created
Definition: carma_context.h:104
carma_magma_syevd
int carma_magma_syevd(char jobz, CarmaObj< T > *mat, CarmaHostObj< T > *eigenvals)
CarmaHostObj
this class provides wrappers to the generic carma host object
Definition: carma_host_obj.h:68
test_rtcFHF.frame
frame
Definition: test_rtcFHF.py:36
correlation_study.k
int k
Definition: correlation_study.py:482
CarmaObj::reset
int reset()
Definition: carma_obj.h:270
layers_test.name
string name
Definition: layers_test.py:21