1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatGetRowIJ for SeqAIJ, SeqBAIJ and SeqSBAIJ\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown PetscErrorCode DumpCSR(Mat A,PetscInt shift,PetscBool symmetric,PetscBool compressed) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown MatType type; 9c4762a1bSJed Brown PetscInt i,j,nr,bs = 1; 10c4762a1bSJed Brown const PetscInt *ia,*ja; 11c4762a1bSJed Brown PetscBool done,isseqbaij,isseqsbaij; 12c4762a1bSJed Brown 13c4762a1bSJed Brown PetscFunctionBeginUser; 14*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQBAIJ,&isseqbaij)); 15*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQSBAIJ,&isseqsbaij)); 16c4762a1bSJed Brown if (isseqbaij || isseqsbaij) { 17*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A,&bs)); 18c4762a1bSJed Brown } 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetType(A,&type)); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(A,shift,symmetric,compressed,&nr,&ia,&ja,&done)); 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"===========================================================\n")); 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"CSR for %s: shift %" PetscInt_FMT " symmetric %" PetscInt_FMT " compressed %" PetscInt_FMT "\n",type,shift,(PetscInt)symmetric,(PetscInt)compressed)); 23c4762a1bSJed Brown for (i=0;i<nr;i++) { 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT ":",i+shift)); 25c4762a1bSJed Brown for (j=ia[i];j<ia[i+1];j++) { 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," %" PetscInt_FMT,ja[j-shift])); 27c4762a1bSJed Brown } 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n")); 29c4762a1bSJed Brown } 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(A,shift,symmetric,compressed,&nr,&ia,&ja,&done)); 31c4762a1bSJed Brown PetscFunctionReturn(0); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34c4762a1bSJed Brown int main(int argc,char **args) 35c4762a1bSJed Brown { 36c4762a1bSJed Brown Mat A,B,C; 37c4762a1bSJed Brown PetscInt i,j,k,m = 3,n = 3,bs = 1; 38c4762a1bSJed Brown PetscErrorCode ierr; 39c4762a1bSJed Brown PetscMPIInt size; 40c4762a1bSJed Brown 41c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 42*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 432c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 47c4762a1bSJed Brown /* adjust sizes by block size */ 48c4762a1bSJed Brown if (m%bs) m += bs-m%bs; 49c4762a1bSJed Brown if (n%bs) n += bs-n%bs; 50c4762a1bSJed Brown 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(A,bs)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSEQAIJ)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&B)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(B,bs)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATSEQBAIJ)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&C)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSize(C,bs)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,MATSEQSBAIJ)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown for (i=0; i<m; i++) { 69c4762a1bSJed Brown for (j=0; j<n; j++) { 70c4762a1bSJed Brown 71c4762a1bSJed Brown PetscScalar v = -1.0; 72c4762a1bSJed Brown PetscInt Ii = j + n*i,J; 73c4762a1bSJed Brown J = Ii - n; 74c4762a1bSJed Brown if (J>=0) { 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown J = Ii + n; 80c4762a1bSJed Brown if (J<m*n) { 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown J = Ii - 1; 86c4762a1bSJed Brown if (J>=0) { 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown J = Ii + 1; 92c4762a1bSJed Brown if (J<m*n) { 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown v = 4.0; 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown } 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* test MatGetRowIJ for the three Mat types */ 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,NULL)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,NULL)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,NULL)); 114c4762a1bSJed Brown for (i=0;i<2;i++) { 115c4762a1bSJed Brown PetscInt shift = i; 116c4762a1bSJed Brown for (j=0;j<2;j++) { 117c4762a1bSJed Brown PetscBool symmetric = ((j>0) ? PETSC_FALSE : PETSC_TRUE); 118c4762a1bSJed Brown for (k=0;k<2;k++) { 119c4762a1bSJed Brown PetscBool compressed = ((k>0) ? PETSC_FALSE : PETSC_TRUE); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(DumpCSR(A,shift,symmetric,compressed)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(DumpCSR(B,shift,symmetric,compressed)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(DumpCSR(C,shift,symmetric,compressed)); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 125c4762a1bSJed Brown } 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 129c4762a1bSJed Brown ierr = PetscFinalize(); 130c4762a1bSJed Brown return ierr; 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown /*TEST 134c4762a1bSJed Brown 135c4762a1bSJed Brown test: 136c4762a1bSJed Brown 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown suffix: 2 139c4762a1bSJed Brown args: -bs 2 140c4762a1bSJed Brown 141c4762a1bSJed Brown TEST*/ 142