#ifndef H_DKS_OPAL #define H_DKS_OPAL #include #include "AutoTuning/DKSAutoTuning.h" #include "DKSBase.h" #include "DKSFFT.h" #include "DKSDefinitions.h" #include "Algorithms/GreensFunction.h" #include "Algorithms/CollimatorPhysics.h" #include "Algorithms/FFT.h" #ifdef DKS_AMD #include "OpenCL/OpenCLFFT.h" #include "OpenCL/OpenCLGreensFunction.h" #include "OpenCL/OpenCLCollimatorPhysics.h" #endif #ifdef DKS_CUDA #include "CUDA/CudaFFT.cuh" #include "CUDA/CudaGreensFunction.cuh" #include "CUDA/CudaCollimatorPhysics.cuh" #endif #ifdef DKS_MIC #include "MIC/MICFFT.h" #include "MIC/MICGreensFunction.hpp" #include "MIC/MICCollimatorPhysics.h" #endif /** * API to handle OPAL calls to DKS library. * Gives access to DKSCollimatorPhysics, GreensFunction and DKSFFT, as well as all the DKSBase * functions. */ class DKSOPAL : public DKSFFT { private: DKSCollimatorPhysics *dkscol; GreensFunction *dksgreens; int setupOPAL(); public: DKSOPAL(); DKSOPAL(const char* api_name, const char* device_name); ~DKSOPAL(); int initDevice(); /////////////////////////////////////////////// ///////Function library part of dksbase//////// /////////////////////////////////////////////// /** * Integrated greens function from OPAL FFTPoissonsolver.cpp put on device. * For specifics check OPAL docs. * TODO: opencl and mic implementations. */ int callGreensIntegral(void *tmp_ptr, int I, int J, int K, int NI, int NJ, double hz_m0, double hz_m1, double hz_m2, int streamId = -1); /** * Integrated greens function from OPAL FFTPoissonsolver.cpp put on device. * For specifics check OPAL docs. * TODO: opencl and mic implementations. */ int callGreensIntegration(void *mem_ptr, void *tmp_ptr, int I, int J, int K, int streamId = -1); /** * Integrated greens function from OPAL FFTPoissonsolver.cpp put on device. * For specifics check OPAL docs. * TODO: opencl and mic implementations. */ int callMirrorRhoField(void *mem_ptr, int I, int J, int K, int streamId = -1); /** * Element by element multiplication. * Multiplies each element of mem_ptr1 with corresponding element of mem_ptr2, size specifies * the number of elements in mem_ptr1 and mem_ptr2 to use. Results are put in mem_ptr1. * TODO: opencl and mic implementations. */ int callMultiplyComplexFields(void *mem_ptr1, void *mem_ptr2, int size, int streamId = -1); /** * Monte carlo code for the degrader from OPAL classic/5.0/src/Solvers/CollimatorPhysics.cpp on device. * For specifics check OPAL docs and CudaCollimatorPhysics class documentation. * TODO: opencl and mic implementations. */ int callCollimatorPhysics(void *mem_ptr, void *par_ptr, int numparticles, int numparams, int &numaddback, int &numdead, bool enableRutherfordScattering = true); /** * Monte carlo code for the degrader from OPAL classic/5.0/src/Solvers/CollimatorPhysics.cpp on device. * For specifics check OPAL docs and CudaCollimatorPhysics class documentation. * TODO: opencl and mic implementations. */ int callCollimatorPhysics2(void *mem_ptr, void *par_ptr, int numparticles, bool enableRutherfordScattering = true); /** * Monte carlo code for the degrader from OPAL classic/5.0/src/Solvers/CollimatorPhysics.cpp on device. * For specifics check OPAL docs and CudaCollimatorPhysics class documentation. * Test function for the MIC to test SoA layout vs AoS layout used in previous versions */ int callCollimatorPhysicsSoA(void *label_ptr, void *localID_ptr, void *rx_ptr, void *ry_ptr, void *rz_ptr, void *px_ptr, void *py_ptr, void *pz_ptr, void *par_ptr, int numparticles); /** * Monte carlo code for the degrader from OPAL classic/5.0/src/Solvers/CollimatorPhysics.cpp on device. * For specifics check OPAL docs and CudaCollimatorPhysics class documentation. * TODO: opencl and mic implementations. */ int callCollimatorPhysicsSort(void *mem_ptr, int numparticles, int &numaddback); /** * Monte carlo code for the degrader from OPAL classic/5.0/src/Solvers/CollimatorPhysics.cpp on device. * For specifics check OPAL docs and CudaCollimatorPhysics class documentation. * TODO: opencl and mic implementations. */ int callCollimatorPhysicsSortSoA(void *label_ptr, void *localID_ptr, void *rx_ptr, void *ry_ptr, void *rz_ptr, void *px_ptr, void *py_ptr, void *pz_ptr, void *par_ptr, int numparticles, int &numaddback); /** * Integration code from ParallelTTracker from OPAL. * For specifics check OPAL docs and CudaCollimatorPhysics class docs */ int callParallelTTrackerPush(void *r_ptr, void *p_ptr, int npart, void *dt_ptr, double dt, double c, bool usedt = false, int streamId = -1); /** * Integration code from ParallelTTracker from OPAL. * For specifics check OPAL docs and CudaCollimatorPhysics class docs */ int callParallelTTrackerPushTransform(void *x_ptr, void *p_ptr, void *lastSec_ptr, void *orient_ptr, int npart, int nsec, void *dt_ptr, double dt, double c, bool usedt = false, int streamId = -1); /** * Integration code from ParallelTTracker from OPAL. * For specifics check OPAL docs and CudaCollimatorPhysics class docs */ int callParallelTTrackerPush(void *r_ptr, void *p_ptr, void *dt_ptr, int npart, double c, int streamId = -1); /** * Integration code from ParallelTTracker from OPAL. * For specifics check OPAL docs and CudaCollimatorPhysics class docs */ int callParallelTTrackerKick(void *r_ptr, void *p_ptr, void *ef_ptr, void *bf_ptr, void *dt_ptr, double charge, double mass, int npart, double c, int streamId = -1); }; #endif