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