xref: /petsc/src/mat/tests/ex26.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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;
149566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQBAIJ,&isseqbaij));
159566063dSJacob Faibussowitsch   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQSBAIJ,&isseqsbaij));
16c4762a1bSJed Brown   if (isseqbaij || isseqsbaij) {
179566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(A,&bs));
18c4762a1bSJed Brown   }
199566063dSJacob Faibussowitsch   PetscCall(MatGetType(A,&type));
209566063dSJacob Faibussowitsch   PetscCall(MatGetRowIJ(A,shift,symmetric,compressed,&nr,&ia,&ja,&done));
219566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"===========================================================\n"));
229566063dSJacob 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++) {
249566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT ":",i+shift));
25c4762a1bSJed Brown     for (j=ia[i];j<ia[i+1];j++) {
269566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF," %" PetscInt_FMT,ja[j-shift]));
27c4762a1bSJed Brown     }
289566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n"));
29c4762a1bSJed Brown   }
309566063dSJacob 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*327415f7SBarry Smith   PetscFunctionBeginUser;
419566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
43be096a46SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
47c4762a1bSJed Brown   /* adjust sizes by block size */
48c4762a1bSJed Brown   if (m%bs) m += bs-m%bs;
49c4762a1bSJed Brown   if (n%bs) n += bs-n%bs;
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&A));
529566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
539566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(A,bs));
549566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATSEQAIJ));
559566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
569566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&B));
579566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
589566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(B,bs));
599566063dSJacob Faibussowitsch   PetscCall(MatSetType(B,MATSEQBAIJ));
609566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
619566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&C));
629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
639566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(C,bs));
649566063dSJacob Faibussowitsch   PetscCall(MatSetType(C,MATSEQSBAIJ));
659566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
669566063dSJacob Faibussowitsch   PetscCall(MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   for (i=0; i<m; i++) {
69c4762a1bSJed Brown     for (j=0; j<n; j++) {
70c4762a1bSJed Brown 
71c4762a1bSJed Brown       PetscScalar v = -1.0;
72c4762a1bSJed Brown       PetscInt    Ii = j + n*i,J;
73c4762a1bSJed Brown       J = Ii - n;
74c4762a1bSJed Brown       if (J>=0)  {
759566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
769566063dSJacob Faibussowitsch         PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
779566063dSJacob Faibussowitsch         PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
78c4762a1bSJed Brown       }
79c4762a1bSJed Brown       J = Ii + n;
80c4762a1bSJed Brown       if (J<m*n) {
819566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
829566063dSJacob Faibussowitsch         PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
839566063dSJacob Faibussowitsch         PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
84c4762a1bSJed Brown       }
85c4762a1bSJed Brown       J = Ii - 1;
86c4762a1bSJed Brown       if (J>=0)  {
879566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
889566063dSJacob Faibussowitsch         PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
899566063dSJacob Faibussowitsch         PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
90c4762a1bSJed Brown       }
91c4762a1bSJed Brown       J = Ii + 1;
92c4762a1bSJed Brown       if (J<m*n) {
939566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
949566063dSJacob Faibussowitsch         PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
959566063dSJacob Faibussowitsch         PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
96c4762a1bSJed Brown       }
97c4762a1bSJed Brown       v = 4.0;
989566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
999566063dSJacob Faibussowitsch       PetscCall(MatSetValues(B,1,&Ii,1,&Ii,&v,INSERT_VALUES));
1009566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
101c4762a1bSJed Brown     }
102c4762a1bSJed Brown   }
1039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
1059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
1089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* test MatGetRowIJ for the three Mat types */
1119566063dSJacob Faibussowitsch   PetscCall(MatView(A,NULL));
1129566063dSJacob Faibussowitsch   PetscCall(MatView(B,NULL));
1139566063dSJacob Faibussowitsch   PetscCall(MatView(C,NULL));
114c4762a1bSJed Brown   for (i=0;i<2;i++) {
115c4762a1bSJed Brown     PetscInt shift = i;
116c4762a1bSJed Brown     for (j=0;j<2;j++) {
117c4762a1bSJed Brown       PetscBool symmetric = ((j>0) ? PETSC_FALSE : PETSC_TRUE);
118c4762a1bSJed Brown       for (k=0;k<2;k++) {
119c4762a1bSJed Brown         PetscBool compressed = ((k>0) ? PETSC_FALSE : PETSC_TRUE);
1209566063dSJacob Faibussowitsch         PetscCall(DumpCSR(A,shift,symmetric,compressed));
1219566063dSJacob Faibussowitsch         PetscCall(DumpCSR(B,shift,symmetric,compressed));
1229566063dSJacob Faibussowitsch         PetscCall(DumpCSR(C,shift,symmetric,compressed));
123c4762a1bSJed Brown       }
124c4762a1bSJed Brown     }
125c4762a1bSJed Brown   }
1269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1299566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
130b122ec5aSJacob Faibussowitsch   return 0;
131c4762a1bSJed Brown }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown /*TEST
134c4762a1bSJed Brown 
135c4762a1bSJed Brown    test:
136c4762a1bSJed Brown 
137c4762a1bSJed Brown    test:
138c4762a1bSJed Brown       suffix: 2
139c4762a1bSJed Brown       args: -bs 2
140c4762a1bSJed Brown 
141c4762a1bSJed Brown TEST*/
142