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