xref: /petsc/src/mat/tests/ex26.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQBAIJ,&isseqbaij));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQSBAIJ,&isseqsbaij));
16c4762a1bSJed Brown   if (isseqbaij || isseqsbaij) {
175f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetBlockSize(A,&bs));
18c4762a1bSJed Brown   }
195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetType(A,&type));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetRowIJ(A,shift,symmetric,compressed,&nr,&ia,&ja,&done));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"===========================================================\n"));
225f80ce2aSJacob 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++) {
245f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT ":",i+shift));
25c4762a1bSJed Brown     for (j=ia[i];j<ia[i+1];j++) {
265f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF," %" PetscInt_FMT,ja[j-shift]));
27c4762a1bSJed Brown     }
285f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
29c4762a1bSJed Brown   }
305f80ce2aSJacob 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   PetscMPIInt    size;
39c4762a1bSJed Brown 
40*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
415f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
422c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(A,bs));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATSEQAIJ));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&B));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(B,bs));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATSEQBAIJ));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(B));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(C,bs));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(C,MATSEQSBAIJ));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(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)  {
745f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
755f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
765f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
77c4762a1bSJed Brown       }
78c4762a1bSJed Brown       J = Ii + n;
79c4762a1bSJed Brown       if (J<m*n) {
805f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
815f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
825f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
83c4762a1bSJed Brown       }
84c4762a1bSJed Brown       J = Ii - 1;
85c4762a1bSJed Brown       if (J>=0)  {
865f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
875f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
885f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
89c4762a1bSJed Brown       }
90c4762a1bSJed Brown       J = Ii + 1;
91c4762a1bSJed Brown       if (J<m*n) {
925f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
935f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
945f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
95c4762a1bSJed Brown       }
96c4762a1bSJed Brown       v = 4.0;
975f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
985f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(B,1,&Ii,1,&Ii,&v,INSERT_VALUES));
995f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
100c4762a1bSJed Brown     }
101c4762a1bSJed Brown   }
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* test MatGetRowIJ for the three Mat types */
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,NULL));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,NULL));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(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);
1195f80ce2aSJacob Faibussowitsch         CHKERRQ(DumpCSR(A,shift,symmetric,compressed));
1205f80ce2aSJacob Faibussowitsch         CHKERRQ(DumpCSR(B,shift,symmetric,compressed));
1215f80ce2aSJacob Faibussowitsch         CHKERRQ(DumpCSR(C,shift,symmetric,compressed));
122c4762a1bSJed Brown       }
123c4762a1bSJed Brown     }
124c4762a1bSJed Brown   }
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
128*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
129*b122ec5aSJacob 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