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