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