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