1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MATSEQDENSECUDA\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **argv) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A,AC,B; 9c4762a1bSJed Brown PetscInt m = 10,n = 10; 10c4762a1bSJed Brown PetscReal r,tol = 10*PETSC_SMALL; 11c4762a1bSJed Brown 12*327415f7SBarry Smith PetscFunctionBeginUser; 139566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*) 0,help)); 149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 159566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 169566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n)); 189566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQDENSE)); 199566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 209566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(A,NULL)); 219566063dSJacob Faibussowitsch PetscCall(MatSetRandom(A,NULL)); 22c20d7725SJed Brown #if 0 23c20d7725SJed Brown PetscInt i,j; 24c20d7725SJed Brown PetscScalar val; 25c20d7725SJed Brown for (i=0; i<m; i++) { 26c20d7725SJed Brown for (j=0; j<n; j++) { 27c20d7725SJed Brown val = (PetscScalar)(i+j); 289566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&val,INSERT_VALUES)); 29c20d7725SJed Brown } 30c20d7725SJed Brown } 319566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 329566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 33c20d7725SJed Brown #endif 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* Create a CUDA version of A */ 36c20d7725SJed Brown #if defined(PETSC_HAVE_CUDA) 379566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSEQDENSECUDA,MAT_INITIAL_MATRIX,&AC)); 38c20d7725SJed Brown #else 399566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&AC)); 40c20d7725SJed Brown #endif 419566063dSJacob Faibussowitsch PetscCall(MatDuplicate(AC,MAT_COPY_VALUES,&B)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* full CUDA AXPY */ 449566063dSJacob Faibussowitsch PetscCall(MatAXPY(B,-1.0,AC,SAME_NONZERO_PATTERN)); 459566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_INFINITY,&r)); 4608401ef6SPierre Jolivet PetscCheck(r == 0.0,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatDuplicate + MatCopy + MatAXPY %g",(double)r); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* test Copy */ 499566063dSJacob Faibussowitsch PetscCall(MatCopy(AC,B,SAME_NONZERO_PATTERN)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* call MatAXPY_Basic since B is CUDA, A is CPU, */ 529566063dSJacob Faibussowitsch PetscCall(MatAXPY(B,-1.0,A,SAME_NONZERO_PATTERN)); 539566063dSJacob Faibussowitsch PetscCall(MatNorm(B,NORM_INFINITY,&r)); 5408401ef6SPierre Jolivet PetscCheck(r == 0.0,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatDuplicate + MatCopy + MatAXPY_Basic %g",(double)r); 55c4762a1bSJed Brown 56c4762a1bSJed Brown if (m == n) { 57c4762a1bSJed Brown Mat B1,B2; 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(MatCopy(AC,B,SAME_NONZERO_PATTERN)); 60c4762a1bSJed Brown /* full CUDA PtAP */ 619566063dSJacob Faibussowitsch PetscCall(MatPtAP(B,AC,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B1)); 62c20d7725SJed Brown 63c4762a1bSJed Brown /* CPU PtAP since A is on the CPU only */ 649566063dSJacob Faibussowitsch PetscCall(MatPtAP(B,A,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B2)); 65c20d7725SJed Brown 669566063dSJacob Faibussowitsch PetscCall(MatAXPY(B2,-1.0,B1,SAME_NONZERO_PATTERN)); 679566063dSJacob Faibussowitsch PetscCall(MatNorm(B2,NORM_INFINITY,&r)); 6808401ef6SPierre Jolivet PetscCheck(r <= tol,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatPtAP %g",(double)r); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* test reuse */ 719566063dSJacob Faibussowitsch PetscCall(MatPtAP(B,AC,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B1)); 729566063dSJacob Faibussowitsch PetscCall(MatPtAP(B,A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B2)); 739566063dSJacob Faibussowitsch PetscCall(MatAXPY(B2,-1.0,B1,SAME_NONZERO_PATTERN)); 749566063dSJacob Faibussowitsch PetscCall(MatNorm(B2,NORM_INFINITY,&r)); 7508401ef6SPierre Jolivet PetscCheck(r <= tol,PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"Error MatPtAP %g",(double)r); 76c4762a1bSJed Brown 779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B1)); 789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AC)); 839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 849566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 85b122ec5aSJacob Faibussowitsch return 0; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown /*TEST 89c4762a1bSJed Brown 90c4762a1bSJed Brown build: 91c4762a1bSJed Brown requires: cuda 92c4762a1bSJed Brown 93c4762a1bSJed Brown test: 94c4762a1bSJed Brown output_file: output/ex32_1.out 95c4762a1bSJed Brown args: -m {{3 5 12}} -n {{3 5 12}} 96c4762a1bSJed Brown suffix: seqdensecuda 97c4762a1bSJed Brown 98c4762a1bSJed Brown TEST*/ 99