Actual source code: cusparsematimpl.h

  1: #if !defined(CUSPARSEMATIMPL)
  2: #define CUSPARSEMATIMPL

  4: #include <petscpkg_version.h>
  5: #include <petsc/private/cudavecimpl.h>
  6: #include <petscaijdevice.h>

  8: #include <cusparse_v2.h>

 10: #include <algorithm>
 11: #include <vector>

 13: #include <thrust/device_vector.h>
 14: #include <thrust/device_ptr.h>
 15: #include <thrust/device_malloc_allocator.h>
 16: #include <thrust/transform.h>
 17: #include <thrust/functional.h>
 18: #include <thrust/sequence.h>
 19: #include <thrust/system/system_error.h>

 21: #define PetscStackCallThrust(body) do {                                     \
 22:     try {                                                                   \
 23:       body;                                                                 \
 24:     } catch(thrust::system_error& e) {                                      \
 25:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Thrust %s",e.what());\
 26:     }                                                                       \
 27:   } while (0)

 29: #if defined(PETSC_USE_COMPLEX)
 30:   #if defined(PETSC_USE_REAL_SINGLE)
 31:     const cuComplex PETSC_CUSPARSE_ONE        = {1.0f, 0.0f};
 32:     const cuComplex PETSC_CUSPARSE_ZERO       = {0.0f, 0.0f};
 33:   #elif defined(PETSC_USE_REAL_DOUBLE)
 34:     const cuDoubleComplex PETSC_CUSPARSE_ONE  = {1.0, 0.0};
 35:     const cuDoubleComplex PETSC_CUSPARSE_ZERO = {0.0, 0.0};
 36:   #endif
 37: #else
 38:   const PetscScalar PETSC_CUSPARSE_ONE        = 1.0;
 39:   const PetscScalar PETSC_CUSPARSE_ZERO       = 0.0;
 40: #endif

 42: #if PETSC_PKG_CUDA_VERSION_GE(9,0,0)
 43:   #define cusparse_create_analysis_info  cusparseCreateCsrsv2Info
 44:   #define cusparse_destroy_analysis_info cusparseDestroyCsrsv2Info
 45:   #if defined(PETSC_USE_COMPLEX)
 46:     #if defined(PETSC_USE_REAL_SINGLE)
 47:       #define cusparse_get_svbuffsize(a,b,c,d,e,f,g,h,i,j) cusparseCcsrsv2_bufferSize(a,b,c,d,e,(cuComplex*)(f),g,h,i,j)
 48:       #define cusparse_analysis(a,b,c,d,e,f,g,h,i,j,k)     cusparseCcsrsv2_analysis(a,b,c,d,e,(const cuComplex*)(f),g,h,i,j,k)
 49:       #define cusparse_solve(a,b,c,d,e,f,g,h,i,j,k,l,m,n)  cusparseCcsrsv2_solve(a,b,c,d,(const cuComplex*)(e),f,(const cuComplex*)(g),h,i,j,(const cuComplex*)(k),(cuComplex*)(l),m,n)
 50:     #elif defined(PETSC_USE_REAL_DOUBLE)
 51:       #define cusparse_get_svbuffsize(a,b,c,d,e,f,g,h,i,j) cusparseZcsrsv2_bufferSize(a,b,c,d,e,(cuDoubleComplex*)(f),g,h,i,j)
 52:       #define cusparse_analysis(a,b,c,d,e,f,g,h,i,j,k)     cusparseZcsrsv2_analysis(a,b,c,d,e,(const cuDoubleComplex*)(f),g,h,i,j,k)
 53:       #define cusparse_solve(a,b,c,d,e,f,g,h,i,j,k,l,m,n)  cusparseZcsrsv2_solve(a,b,c,d,(const cuDoubleComplex*)(e),f,(const cuDoubleComplex*)(g),h,i,j,(const cuDoubleComplex*)(k),(cuDoubleComplex*)(l),m,n)
 54:     #endif
 55:   #else /* not complex */
 56:     #if defined(PETSC_USE_REAL_SINGLE)
 57:       #define cusparse_get_svbuffsize cusparseScsrsv2_bufferSize
 58:       #define cusparse_analysis       cusparseScsrsv2_analysis
 59:       #define cusparse_solve          cusparseScsrsv2_solve
 60:     #elif defined(PETSC_USE_REAL_DOUBLE)
 61:       #define cusparse_get_svbuffsize cusparseDcsrsv2_bufferSize
 62:       #define cusparse_analysis       cusparseDcsrsv2_analysis
 63:       #define cusparse_solve          cusparseDcsrsv2_solve
 64:     #endif
 65:   #endif
 66: #else
 67:   #define cusparse_create_analysis_info  cusparseCreateSolveAnalysisInfo
 68:   #define cusparse_destroy_analysis_info cusparseDestroySolveAnalysisInfo
 69:   #if defined(PETSC_USE_COMPLEX)
 70:     #if defined(PETSC_USE_REAL_SINGLE)
 71:       #define cusparse_solve(a,b,c,d,e,f,g,h,i,j,k) cusparseCcsrsv_solve((a),(b),(c),(cuComplex*)(d),(e),(cuComplex*)(f),(g),(h),(i),(cuComplex*)(j),(cuComplex*)(k))
 72:       #define cusparse_analysis(a,b,c,d,e,f,g,h,i)  cusparseCcsrsv_analysis((a),(b),(c),(d),(e),(cuComplex*)(f),(g),(h),(i))
 73:     #elif defined(PETSC_USE_REAL_DOUBLE)
 74:       #define cusparse_solve(a,b,c,d,e,f,g,h,i,j,k) cusparseZcsrsv_solve((a),(b),(c),(cuDoubleComplex*)(d),(e),(cuDoubleComplex*)(f),(g),(h),(i),(cuDoubleComplex*)(j),(cuDoubleComplex*)(k))
 75:       #define cusparse_analysis(a,b,c,d,e,f,g,h,i)  cusparseZcsrsv_analysis((a),(b),(c),(d),(e),(cuDoubleComplex*)(f),(g),(h),(i))
 76:     #endif
 77:   #else /* not complex */
 78:     #if defined(PETSC_USE_REAL_SINGLE)
 79:       #define cusparse_solve    cusparseScsrsv_solve
 80:       #define cusparse_analysis cusparseScsrsv_analysis
 81:     #elif defined(PETSC_USE_REAL_DOUBLE)
 82:       #define cusparse_solve    cusparseDcsrsv_solve
 83:       #define cusparse_analysis cusparseDcsrsv_analysis
 84:     #endif
 85:   #endif
 86: #endif

 88: #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
 89:   #define cusparse_csr2csc cusparseCsr2cscEx2
 90:   #if defined(PETSC_USE_COMPLEX)
 91:     #if defined(PETSC_USE_REAL_SINGLE)
 92:       #define cusparse_scalartype CUDA_C_32F
 93:       #define cusparse_csr_spgeam(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t)            cusparseCcsrgeam2(a,b,c,(cuComplex*)d,e,f,(cuComplex*)g,h,i,(cuComplex*)j,k,l,(cuComplex*)m,n,o,p,(cuComplex*)q,r,s,t)
 94:       #define cusparse_csr_spgeam_bufferSize(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t) cusparseCcsrgeam2_bufferSizeExt(a,b,c,(cuComplex*)d,e,f,(cuComplex*)g,h,i,(cuComplex*)j,k,l,(cuComplex*)m,n,o,p,(cuComplex*)q,r,s,t)
 95:     #elif defined(PETSC_USE_REAL_DOUBLE)
 96:       #define cusparse_scalartype CUDA_C_64F
 97:       #define cusparse_csr_spgeam(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t)            cusparseZcsrgeam2(a,b,c,(cuDoubleComplex*)d,e,f,(cuDoubleComplex*)g,h,i,(cuDoubleComplex*)j,k,l,(cuDoubleComplex*)m,n,o,p,(cuDoubleComplex*)q,r,s,t)
 98:       #define cusparse_csr_spgeam_bufferSize(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t) cusparseZcsrgeam2_bufferSizeExt(a,b,c,(cuDoubleComplex*)d,e,f,(cuDoubleComplex*)g,h,i,(cuDoubleComplex*)j,k,l,(cuDoubleComplex*)m,n,o,p,(cuDoubleComplex*)q,r,s,t)
 99:     #endif
100:   #else /* not complex */
101:     #if defined(PETSC_USE_REAL_SINGLE)
102:       #define cusparse_scalartype CUDA_R_32F
103:       #define cusparse_csr_spgeam            cusparseScsrgeam2
104:       #define cusparse_csr_spgeam_bufferSize cusparseScsrgeam2_bufferSizeExt
105:     #elif defined(PETSC_USE_REAL_DOUBLE)
106:       #define cusparse_scalartype CUDA_R_64F
107:       #define cusparse_csr_spgeam            cusparseDcsrgeam2
108:       #define cusparse_csr_spgeam_bufferSize cusparseDcsrgeam2_bufferSizeExt
109:     #endif
110:   #endif
111: #else
112:   #if defined(PETSC_USE_COMPLEX)
113:     #if defined(PETSC_USE_REAL_SINGLE)
114:       #define cusparse_csr_spmv(a,b,c,d,e,f,g,h,i,j,k,l,m)       cusparseCcsrmv((a),(b),(c),(d),(e),(cuComplex*)(f),(g),(cuComplex*)(h),(i),(j),(cuComplex*)(k),(cuComplex*)(l),(cuComplex*)(m))
115:       #define cusparse_csr_spmm(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p) cusparseCcsrmm((a),(b),(c),(d),(e),(f),(cuComplex*)(g),(h),(cuComplex*)(i),(j),(k),(cuComplex*)(l),(m),(cuComplex*)(n),(cuComplex*)(o),(p))
116:       #define cusparse_csr2csc(a,b,c,d,e,f,g,h,i,j,k,l)          cusparseCcsr2csc((a),(b),(c),(d),(cuComplex*)(e),(f),(g),(cuComplex*)(h),(i),(j),(k),(l))
117:       #define cusparse_hyb_spmv(a,b,c,d,e,f,g,h)                 cusparseChybmv((a),(b),(cuComplex*)(c),(d),(e),(cuComplex*)(f),(cuComplex*)(g),(cuComplex*)(h))
118:       #define cusparse_csr2hyb(a,b,c,d,e,f,g,h,i,j)              cusparseCcsr2hyb((a),(b),(c),(d),(cuComplex*)(e),(f),(g),(h),(i),(j))
119:       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseChyb2csr((a),(b),(c),(cuComplex*)(d),(e),(f))
120:       #define cusparse_csr_spgemm(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t) cusparseCcsrgemm(a,b,c,d,e,f,g,h,(cuComplex*)i,j,k,l,m,(cuComplex*)n,o,p,q,(cuComplex*)r,s,t)
121:       #define cusparse_csr_spgeam(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s)   cusparseCcsrgeam(a,b,c,(cuComplex*)d,e,f,(cuComplex*)g,h,i,(cuComplex*)j,k,l,(cuComplex*)m,n,o,p,(cuComplex*)q,r,s)
122:     #elif defined(PETSC_USE_REAL_DOUBLE)
123:       #define cusparse_csr_spmv(a,b,c,d,e,f,g,h,i,j,k,l,m)       cusparseZcsrmv((a),(b),(c),(d),(e),(cuDoubleComplex*)(f),(g),(cuDoubleComplex*)(h),(i),(j),(cuDoubleComplex*)(k),(cuDoubleComplex*)(l),(cuDoubleComplex*)(m))
124:       #define cusparse_csr_spmm(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p) cusparseZcsrmm((a),(b),(c),(d),(e),(f),(cuDoubleComplex*)(g),(h),(cuDoubleComplex*)(i),(j),(k),(cuDoubleComplex*)(l),(m),(cuDoubleComplex*)(n),(cuDoubleComplex*)(o),(p))
125:       #define cusparse_csr2csc(a,b,c,d,e,f,g,h,i,j,k,l)          cusparseZcsr2csc((a),(b),(c),(d),(cuDoubleComplex*)(e),(f),(g),(cuDoubleComplex*)(h),(i),(j),(k),(l))
126:       #define cusparse_hyb_spmv(a,b,c,d,e,f,g,h)                 cusparseZhybmv((a),(b),(cuDoubleComplex*)(c),(d),(e),(cuDoubleComplex*)(f),(cuDoubleComplex*)(g),(cuDoubleComplex*)(h))
127:       #define cusparse_csr2hyb(a,b,c,d,e,f,g,h,i,j)              cusparseZcsr2hyb((a),(b),(c),(d),(cuDoubleComplex*)(e),(f),(g),(h),(i),(j))
128:       #define cusparse_hyb2csr(a,b,c,d,e,f)                      cusparseZhyb2csr((a),(b),(c),(cuDoubleComplex*)(d),(e),(f))
129:       #define cusparse_csr_spgemm(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s,t) cusparseZcsrgemm(a,b,c,d,e,f,g,h,(cuDoubleComplex*)i,j,k,l,m,(cuDoubleComplex*)n,o,p,q,(cuDoubleComplex*)r,s,t)
130:       #define cusparse_csr_spgeam(a,b,c,d,e,f,g,h,i,j,k,l,m,n,o,p,q,r,s)   cusparseZcsrgeam(a,b,c,(cuDoubleComplex*)d,e,f,(cuDoubleComplex*)g,h,i,(cuDoubleComplex*)j,k,l,(cuDoubleComplex*)m,n,o,p,(cuDoubleComplex*)q,r,s)
131:     #endif
132:   #else
133:     #if defined(PETSC_USE_REAL_SINGLE)
134:       #define cusparse_csr_spmv cusparseScsrmv
135:       #define cusparse_csr_spmm cusparseScsrmm
136:       #define cusparse_csr2csc  cusparseScsr2csc
137:       #define cusparse_hyb_spmv cusparseShybmv
138:       #define cusparse_csr2hyb  cusparseScsr2hyb
139:       #define cusparse_hyb2csr  cusparseShyb2csr
140:       #define cusparse_csr_spgemm cusparseScsrgemm
141:       #define cusparse_csr_spgeam cusparseScsrgeam
142:     #elif defined(PETSC_USE_REAL_DOUBLE)
143:       #define cusparse_csr_spmv cusparseDcsrmv
144:       #define cusparse_csr_spmm cusparseDcsrmm
145:       #define cusparse_csr2csc  cusparseDcsr2csc
146:       #define cusparse_hyb_spmv cusparseDhybmv
147:       #define cusparse_csr2hyb  cusparseDcsr2hyb
148:       #define cusparse_hyb2csr  cusparseDhyb2csr
149:       #define cusparse_csr_spgemm cusparseDcsrgemm
150:       #define cusparse_csr_spgeam cusparseDcsrgeam
151:     #endif
152:   #endif
153: #endif

155: #define THRUSTINTARRAY32 thrust::device_vector<int>
156: #define THRUSTINTARRAY thrust::device_vector<PetscInt>
157: #define THRUSTARRAY thrust::device_vector<PetscScalar>

159: /* A CSR matrix structure */
160: struct CsrMatrix {
161:   PetscInt         num_rows;
162:   PetscInt         num_cols;
163:   PetscInt         num_entries;
164:   THRUSTINTARRAY32 *row_offsets;
165:   THRUSTINTARRAY32 *column_indices;
166:   THRUSTARRAY      *values;
167: };

169: /* This is struct holding the relevant data needed to a MatSolve */
170: struct Mat_SeqAIJCUSPARSETriFactorStruct {
171:   /* Data needed for triangular solve */
172:   cusparseMatDescr_t          descr;
173:   cusparseOperation_t         solveOp;
174:   CsrMatrix                   *csrMat;
175:  #if PETSC_PKG_CUDA_VERSION_GE(9,0,0)
176:   csrsv2Info_t                solveInfo;
177:  #else
178:   cusparseSolveAnalysisInfo_t solveInfo;
179:  #endif
180:   cusparseSolvePolicy_t       solvePolicy;     /* whether level information is generated and used */
181:   int                         solveBufferSize;
182:   void                        *solveBuffer;
183:   size_t                      csr2cscBufferSize; /* to transpose the triangular factor (only used for CUDA >= 11.0) */
184:   void                        *csr2cscBuffer;
185:   PetscScalar                 *AA_h; /* managed host buffer for moving values to the GPU */
186: };

188: /* This is a larger struct holding all the triangular factors for a solve, transpose solve, and any indices used in a reordering */
189: struct Mat_SeqAIJCUSPARSETriFactors {
190:   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtr; /* pointer for lower triangular (factored matrix) on GPU */
191:   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtr; /* pointer for upper triangular (factored matrix) on GPU */
192:   Mat_SeqAIJCUSPARSETriFactorStruct *loTriFactorPtrTranspose; /* pointer for lower triangular (factored matrix) on GPU for the transpose (useful for BiCG) */
193:   Mat_SeqAIJCUSPARSETriFactorStruct *upTriFactorPtrTranspose; /* pointer for upper triangular (factored matrix) on GPU for the transpose (useful for BiCG)*/
194:   THRUSTINTARRAY                    *rpermIndices;  /* indices used for any reordering */
195:   THRUSTINTARRAY                    *cpermIndices;  /* indices used for any reordering */
196:   THRUSTARRAY                       *workVector;
197:   cusparseHandle_t                  handle;   /* a handle to the cusparse library */
198:   PetscInt                          nnz;      /* number of nonzeros ... need this for accurate logging between ICC and ILU */
199:   PetscScalar                       *a_band_d; /* GPU data for banded CSR LU factorization matrix diag(L)=1 */
200:   int                               *i_band_d; /* this could be optimized away */
201:   cudaDeviceProp                    dev_prop;
202:   PetscBool                         init_dev_prop;
203: };

205: struct Mat_CusparseSpMV {
206:   PetscBool             initialized;    /* Don't rely on spmvBuffer != NULL to test if the struct is initialized, */
207:   size_t                spmvBufferSize; /* since I'm not sure if smvBuffer can be NULL even after cusparseSpMV_bufferSize() */
208:   void                  *spmvBuffer;
209:  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)  /* these are present from CUDA 10.1, but PETSc code makes use of them from CUDA 11 on */
210:   cusparseDnVecDescr_t  vecXDescr,vecYDescr; /* descriptor for the dense vectors in y=op(A)x */
211:  #endif
212: };

214: /* This is struct holding the relevant data needed to a MatMult */
215: struct Mat_SeqAIJCUSPARSEMultStruct {
216:   void               *mat;  /* opaque pointer to a matrix. This could be either a cusparseHybMat_t or a CsrMatrix */
217:   cusparseMatDescr_t descr; /* Data needed to describe the matrix for a multiply */
218:   THRUSTINTARRAY     *cprowIndices;   /* compressed row indices used in the parallel SpMV */
219:   PetscScalar        *alpha_one; /* pointer to a device "scalar" storing the alpha parameter in the SpMV */
220:   PetscScalar        *beta_zero; /* pointer to a device "scalar" storing the beta parameter in the SpMV as zero*/
221:   PetscScalar        *beta_one; /* pointer to a device "scalar" storing the beta parameter in the SpMV as one */
222:  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
223:   cusparseSpMatDescr_t  matDescr;  /* descriptor for the matrix, used by SpMV and SpMM */
224:   Mat_CusparseSpMV      cuSpMV[3]; /* different Mat_CusparseSpMV structs for non-transpose, transpose, conj-transpose */
225:   Mat_SeqAIJCUSPARSEMultStruct() : matDescr(NULL) {
226:     for (int i=0; i<3; i++) cuSpMV[i].initialized = PETSC_FALSE;
227:   }
228:  #endif
229: };

231: /* This is a larger struct holding all the matrices for a SpMV, and SpMV Transpose */
232: struct Mat_SeqAIJCUSPARSE {
233:   Mat_SeqAIJCUSPARSEMultStruct *mat;            /* pointer to the matrix on the GPU */
234:   Mat_SeqAIJCUSPARSEMultStruct *matTranspose;   /* pointer to the matrix on the GPU (for the transpose ... useful for BiCG) */
235:   THRUSTARRAY                  *workVector;     /* pointer to a workvector to which we can copy the relevant indices of a vector we want to multiply */
236:   THRUSTINTARRAY32             *rowoffsets_gpu; /* rowoffsets on GPU in non-compressed-row format. It is used to convert CSR to CSC */
237:   PetscInt                     nrows;           /* number of rows of the matrix seen by GPU */
238:   MatCUSPARSEStorageFormat     format;          /* the storage format for the matrix on the device */
239:   PetscBool                    use_cpu_solve;   /* Use AIJ_Seq (I)LU solve */
240:   cudaStream_t                 stream;          /* a stream for the parallel SpMV ... this is not owned and should not be deleted */
241:   cusparseHandle_t             handle;          /* a handle to the cusparse library ... this may not be owned (if we're working in parallel i.e. multiGPUs) */
242:   PetscObjectState             nonzerostate;    /* track nonzero state to possibly recreate the GPU matrix */
243:  #if PETSC_PKG_CUDA_VERSION_GE(11,0,0)
244:   size_t                       csr2cscBufferSize; /* stuff used to compute the matTranspose above */
245:   void                         *csr2cscBuffer;    /* This is used as a C struct and is calloc'ed by PetscNewLog() */
246:   cusparseCsr2CscAlg_t         csr2cscAlg;        /* algorithms can be selected from command line options */
247:   cusparseSpMVAlg_t            spmvAlg;
248:   cusparseSpMMAlg_t            spmmAlg;
249:  #endif
250:   THRUSTINTARRAY               *csr2csc_i;
251:   PetscSplitCSRDataStructure   deviceMat;       /* Matrix on device for, eg, assembly */

253:   /* Stuff for basic COO support */
254:   THRUSTINTARRAY               *cooPerm;        /* permutation array that sorts the input coo entris by row and col */
255:   THRUSTINTARRAY               *cooPerm_a;      /* ordered array that indicate i-th nonzero (after sorting) is the j-th unique nonzero */

257:   /* Stuff for extended COO support */
258:   PetscBool                    use_extended_coo; /* Use extended COO format */
259:   PetscCount                   *jmap_d; /* perm[disp+jmap[i]..disp+jmap[i+1]) gives indices of entries in v[] associated with i-th nonzero of the matrix */
260:   PetscCount                   *perm_d;

262:   Mat_SeqAIJCUSPARSE() : use_extended_coo(PETSC_FALSE), perm_d(NULL), jmap_d(NULL) {}
263: };

265: PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSECopyToGPU(Mat);
266: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_SeqAIJCUSPARSE_Basic(Mat,PetscCount,const PetscInt[],const PetscInt[]);
267: PETSC_INTERN PetscErrorCode MatSetValuesCOO_SeqAIJCUSPARSE_Basic(Mat,const PetscScalar[],InsertMode);
268: PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSEMergeMats(Mat,Mat,MatReuse,Mat*);
269: PETSC_INTERN PetscErrorCode MatSeqAIJCUSPARSETriFactors_Reset(Mat_SeqAIJCUSPARSETriFactors_p*);

271: static inline bool isCudaMem(const void *data)
272: {
273:   cudaError_t                  cerr;
274:   struct cudaPointerAttributes attr;
275:   enum cudaMemoryType          mtype;
276:   cerr = cudaPointerGetAttributes(&attr,data); /* Do not check error since before CUDA 11.0, passing a host pointer returns cudaErrorInvalidValue */
277:   cudaGetLastError(); /* Reset the last error */
278:   #if (CUDART_VERSION < 10000)
279:     mtype = attr.memoryType;
280:   #else
281:     mtype = attr.type;
282:   #endif
283:   if (cerr == cudaSuccess && mtype == cudaMemoryTypeDevice) return true;
284:   else return false;
285: }

287: #endif