1 2 static char help[] = "Tests for bugs in A->offloadmask consistency for GPU matrices\n\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 Mat A; 9 PetscInt i,j,rstart,rend,m = 3; 10 PetscScalar one = 1.0,zero = 0.0,negativeone = -1.0; 11 PetscReal norm; 12 Vec x,y; 13 PetscErrorCode ierr; 14 15 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 16 ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 17 18 for (i=0; i<2; i++) { 19 /* Create the matrix and set it to contain explicit zero entries on the diagonal. */ 20 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 21 ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*m,m*m);CHKERRQ(ierr); 22 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 23 ierr = MatSetUp(A);CHKERRQ(ierr); 24 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 25 ierr = MatCreateVecs(A,&x,&y);CHKERRQ(ierr); 26 ierr = VecSet(x,one);CHKERRQ(ierr); 27 ierr = VecSet(y,zero);CHKERRQ(ierr); 28 ierr = MatDiagonalSet(A,y,INSERT_VALUES);CHKERRQ(ierr); 29 30 /* Now set A to be the identity using various approaches. 31 * Note that there may be other approaches that should be added here. */ 32 switch (i) { 33 case 0: 34 ierr = MatDiagonalSet(A,x,INSERT_VALUES);CHKERRQ(ierr); 35 break; 36 case 1: 37 for (j=rstart; j<rend; j++) { 38 ierr = MatSetValue(A,j,j,one,INSERT_VALUES);CHKERRQ(ierr); 39 } 40 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42 break; 43 case 2: 44 for (j=rstart; j<rend; j++) { 45 ierr = MatSetValuesRow(A,j,&one);CHKERRQ(ierr); 46 } 47 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 48 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 49 default: 50 break; 51 } 52 53 /* Compute y <- A*x and verify that the difference between y and x is negligible, as it should be since A is the identity. */ 54 ierr = MatMult(A,x,y);CHKERRQ(ierr); 55 ierr = VecAXPY(y,negativeone,x);CHKERRQ(ierr); 56 ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr); 57 if (norm > PETSC_SQRT_MACHINE_EPSILON) { 58 ierr = PetscPrintf(PETSC_COMM_WORLD,"Test %d: Norm of error is %g, but should be near 0.\n",i,(double)norm);CHKERRQ(ierr); 59 } 60 61 ierr = MatDestroy(&A);CHKERRQ(ierr); 62 ierr = VecDestroy(&x);CHKERRQ(ierr); 63 ierr = VecDestroy(&y);CHKERRQ(ierr); 64 } 65 66 ierr = PetscFinalize(); 67 return ierr; 68 } 69 70 /*TEST 71 72 test: 73 suffix: aijviennacl_1 74 nsize: 1 75 args: -mat_type aijviennacl 76 requires: viennacl 77 78 test: 79 suffix: aijviennacl_2 80 nsize: 2 81 args: -mat_type aijviennacl 82 requires: viennacl 83 84 test: 85 suffix: aijcusparse_1 86 nsize: 1 87 args: -mat_type aijcusparse 88 requires: cuda 89 90 test: 91 suffix: aijcusparse_2 92 nsize: 2 93 args: -mat_type aijcusparse 94 requires: cuda 95 TEST*/ 96