Commit c3852083 authored by Giovanni La Mura's avatar Giovanni La Mura
Browse files

Define GPU specific function dependenciess for magma_cms()

parent 4e0bdaa8
Loading
Loading
Loading
Loading
+86 −48
Original line number Diff line number Diff line
@@ -24,24 +24,97 @@
#ifndef INCLUDE_MAGMA_CALLS_H_
#define INCLUDE_MAGMA_CALLS_H_

#define MAX_GPU_STACK_DOUBLE 128

#ifdef USE_TARGET_OFFLOAD
#pragma omp begin declare target device_type(any)
/**
 * \brief C++ porting of CGEV. Get weight of T-matrix element (MAGMA GPU version).
 *
 * \param[in] ipamo: `int`
 * \param[in] mu: `int`
 * \param[in] l: `int`
 * \param[in] m: `int`
 * \return result: `double`
 */
double magma_cgev(int ipamo, int mu, int l, int m);
#pragma omp end declare target

/**
 * \brief Initialize a MAGMA compliant T-matrix leaving it on the GPU.
 *
 * \param vec_am: `magmaDoubleComplex *` Vector form of the matrix.
 * \param c1: `ParticleDescriptor *` Pointer to a ParticleDescriptor instance.
 * \param device_id: `int` ID of the device to build the matrix on.
 * \param[in,out] vec_am: `magmaDoubleComplex *` Vector form of the matrix.
 * \param[in] c1: `ParticleDescriptor *` Pointer to a ParticleDescriptor instance.
 * \param[in] device_id: `int` ID of the device to build the matrix on.
 * \return result: `magma_int_t` A return code (MAGMA_SUCCESS, if everything is fine).
 */
magma_int_t magma_cms(magmaDoubleComplex *vec_am, ParticleDescriptor *c1, int device_id);

/**
 * \brief Compute the transfer vector from N2 to N1.
 *
 * This function computes the transfer vector going from N2 to N1, using either
 * Hankel, Bessel or Bessel from origin functions.
 *
 * \tparam IHI: `int` Function mode.
 * \param[in] ipamo: `int`
 * \param[in] nbl: `int` Block identifier.
 * \param[in] l1: `int` First L quantum number.
 * \param[in] m1: `int` First M quantum number.
 * \param[in] l2: `int` Second L quantum number.
 * \param[in] m2: `int` Second M quantum number.
 * \param[in] rxx: `double *` X coordinates of the spheres.
 * \param[in] ryy: `double *` Y coordinates of the spheres.
 * \param[in] rzz: `double *` Z coordinates of the spheres.
 * \param[in] ind3j: `int *` J vector index look-up table.
 * \param[in] v3j0: `double *` J0 vectors.
 * \param[in] vh: `magmaDoubleComplex *`
 * \param[in] vyhj: `magmaDoubleComplex *`
 * \param[in] vj0: `magmaDoubleComplex *`
 * \param[in] vyj0: `magmaDoubleComplex *`
 * \param[in] vj: `magmaDoubleComplex`
 * \param[in] li: `int` Maximum internal expansion order.
 * \param[in] le: `int` Maximum external expansion order.
 * \param[in,out] rac3j: `double[]` J connection vector.
 * \return result: `magmaDoubleComplex` Matrix element.
 */
#pragma omp begin declare target device_type(any)
template<int IHI> magmaDoubleComplex magma_ghit(
  int ipamo, int nbl, int l1, int m1, int l2, int m2,
  double *rxx, double *ryy, double *rzz, int *ind3j,
  double *v3j0, magmaDoubleComplex *vh, magmaDoubleComplex *vyhj,
  magmaDoubleComplex *vj0, magmaDoubleComplex *vyj0,
  magmaDoubleComplex vj, int li, int le, double rac3j[]
);
#pragma omp end declare target

#pragma omp begin declare target device_type(any)
/**
 * \brief Compute the 3j symbol for Clebsch-Gordan coefficients for JJ transitions.
 *
 * This function calculates the 3j(J,J2,J3;-M2-M3,M2,M3) symbol for the Clebsch-Gordan
 * coefficients. See Appendix a.3.1 in Borghese, Denti & Saija (2007). This function
 * is a specialized version to work on GPUs with MAGMA.
 *
 * \param[in] j2: `int` Value of J2.
 * \param[in] j3: `int` Value of J3.
 * \param[in] m2: `int` Value of M2.
 * \param[in] m3: `int` Value of M3.
 * \param[out] rac3j: `double[]` Vector of 3j symbols.
 */
void magma_r3jjr(int j2, int j3, int m2, int m3, double rac3j[]);
#pragma omp end declare target

#endif // USE_TARGET_OFFLOAD

/**
 * \brief Invert a complex matrix with double precision elements.
 *
 * \param mat: Matrix of complex. The matrix to be inverted.
 * \param n: `np_int` The number of rows and columns of the [n x n] matrix.
 * \param jer: `int &` Reference to an integer return flag.
 * \param device_id: `int` ID of the device for matrix inversion offloading.
 * \param rs: `const RuntimeSettings &` Runtime settings instance.
 * \param[in,out] mat: Matrix of complex. The matrix to be inverted.
 * \param[in] n: `np_int` The number of rows and columns of the [n x n] matrix.
 * \param[out] jer: `int &` Reference to an integer return flag.
 * \param[in] device_id: `int` ID of the device for matrix inversion offloading (optional, default 0).
 * \param[in] rs: `const RuntimeSettings &` Runtime settings instance (optional).
 */
void magma_zinvert(
  dcomplex **mat, np_int n, int &jer, int device_id=0,
@@ -65,11 +138,11 @@ void magma_zinvert(
 * the calling function will be responsible of getting the refined matrix from
 * device back to the host.
 *
 * \param rs: `const RuntimeSettings &` Runtime settings instance. [IN]
 * \param a: `magmaDoubleComplex *` Pointer to the first element of the non-inverted matrix on host. [IN]
 * \param m: `const magma_int_t` Number of rows / columns in a. [IN]
 * \param d_a: `magmaDoubleComplex *` Pointer to the matrix on the GPU. [IN/OUT]
 * \param queue: `magma_queue_t` GPU communication queue. [IN]
 * \param[in] rs: `const RuntimeSettings &` Runtime settings instance.
 * \param[in] a: `magmaDoubleComplex *` Pointer to the first element of the non-inverted matrix on host.
 * \param[in] m: `const magma_int_t` Number of rows / columns in a.
 * \param[in,out] d_a: `magmaDoubleComplex *` Pointer to the matrix on the GPU.
 * \param[in] queue: `magma_queue_t` GPU communication queue.
 * \return err: `magma_int_t` An error code (MAGMA_SUCCESS, if everything was fine).
 */
magma_int_t magma_newton(
@@ -77,39 +150,4 @@ magma_int_t magma_newton(
  magmaDoubleComplex* d_a, magma_queue_t queue
);

/* \brief Invert a complex matrix with double precision elements, applying iterative refinement of the solution
 *
 * call magma_zinvert1() to perform the first matrix inversion, then magma_refine() to do the refinement (only if maxrefiters is >0)
 *
 * \param mat: Matrix of complex. The matrix to be inverted.
 * \param n: `np_int` The number of rows and columns of the [n x n] matrix.
 * \param jer: `int &` Reference to an integer return flag.
 * \param maxrefiters: `int` Maximum number of refinement iterations to apply.
 * \param accuracygoal: `double &` Accuracy to achieve in iterative refinement, defined as the module of the maximum difference between the identity matrix and the matrix product of the (approximate) inverse times the original matrix. On return, it contains the actually achieved accuracy.
 * \param refinemode: `int` Flag to control the refinement mode.
 * \param device_id: `int` ID of the device for matrix inversion offloading.
 * \param output_path: `const string &` Path where the output needs to be placed.
 * \param jxi488: `int` Index of the current wavelength calculation.
 */
//void magma_zinvert_and_refine(dcomplex **mat, np_int n, int &jer, int &maxrefiters, double &accuracygoal, int refinemode, int device_id, const string& output_path, int jxi488);

/* \brief Apply iterative refinement of the solution of a matrix inversion.
 *
 * Iteratively compute and apply a correction to the inverse `inva` of the complex
 * matrix `aorig`, for a maximum number of `maxiters` times, or until achieving a
 * maximum residual better than `accuracygoal`.
 *
 * \param aorig: pointer to the first element of the matrix of complex to be inverted.
 * \param inva: pointer to the first element of inverse.
 * \param n: `np_int` The number of rows and columns of the [n x n] matrices.
 * \param jer: `int &` Reference to an integer return flag.
 * \param maxrefiters: `int` Maximum number of refinement iterations to apply.
 * \param accuracygoal: `double` Accuracy to achieve in iterative refinement, defined as the module of the maximum difference between the identity matrix and the matrix product of the (approximate) inverse times the original matrix. On return, it contains the actually achieved accuracy.
 * \param refinemode: `int` Flag for refinement mode selection.
 * \param device_id: `int` ID of the device for matrix inversion offloading.
 * \param output_path: `const string &` Path where the output needs to be placed.
 * \param jxi488: `int` Index of the current wavelength calculation.
 */
// void magma_refine(dcomplex *aorig, dcomplex *inva, np_int n, int &jer, int &maxrefiters, double &accuracygoal, int refinemode, int device_id, const std::string& output_path, int jxi488);

#endif