xref: /petsc/src/mat/tests/ex26.c (revision daa037dfd3c3bec8dc8659548d2b20b07c1dc6de)
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   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
41   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
42   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
43   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
44   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
45   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
46   /* adjust sizes by block size */
47   if (m%bs) m += bs-m%bs;
48   if (n%bs) n += bs-n%bs;
49 
50   PetscCall(MatCreate(PETSC_COMM_SELF,&A));
51   PetscCall(MatSetSizes(A,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
52   PetscCall(MatSetBlockSize(A,bs));
53   PetscCall(MatSetType(A,MATSEQAIJ));
54   PetscCall(MatSetUp(A));
55   PetscCall(MatCreate(PETSC_COMM_SELF,&B));
56   PetscCall(MatSetSizes(B,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
57   PetscCall(MatSetBlockSize(B,bs));
58   PetscCall(MatSetType(B,MATSEQBAIJ));
59   PetscCall(MatSetUp(B));
60   PetscCall(MatCreate(PETSC_COMM_SELF,&C));
61   PetscCall(MatSetSizes(C,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
62   PetscCall(MatSetBlockSize(C,bs));
63   PetscCall(MatSetType(C,MATSEQSBAIJ));
64   PetscCall(MatSetUp(C));
65   PetscCall(MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE));
66 
67   for (i=0; i<m; i++) {
68     for (j=0; j<n; j++) {
69 
70       PetscScalar v = -1.0;
71       PetscInt    Ii = j + n*i,J;
72       J = Ii - n;
73       if (J>=0)  {
74         PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
75         PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
76         PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
77       }
78       J = Ii + n;
79       if (J<m*n) {
80         PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
81         PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
82         PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
83       }
84       J = Ii - 1;
85       if (J>=0)  {
86         PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
87         PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
88         PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
89       }
90       J = Ii + 1;
91       if (J<m*n) {
92         PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
93         PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
94         PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
95       }
96       v = 4.0;
97       PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
98       PetscCall(MatSetValues(B,1,&Ii,1,&Ii,&v,INSERT_VALUES));
99       PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
100     }
101   }
102   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
103   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
104   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
105   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
106   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
107   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
108 
109   /* test MatGetRowIJ for the three Mat types */
110   PetscCall(MatView(A,NULL));
111   PetscCall(MatView(B,NULL));
112   PetscCall(MatView(C,NULL));
113   for (i=0;i<2;i++) {
114     PetscInt shift = i;
115     for (j=0;j<2;j++) {
116       PetscBool symmetric = ((j>0) ? PETSC_FALSE : PETSC_TRUE);
117       for (k=0;k<2;k++) {
118         PetscBool compressed = ((k>0) ? PETSC_FALSE : PETSC_TRUE);
119         PetscCall(DumpCSR(A,shift,symmetric,compressed));
120         PetscCall(DumpCSR(B,shift,symmetric,compressed));
121         PetscCall(DumpCSR(C,shift,symmetric,compressed));
122       }
123     }
124   }
125   PetscCall(MatDestroy(&A));
126   PetscCall(MatDestroy(&B));
127   PetscCall(MatDestroy(&C));
128   PetscCall(PetscFinalize());
129   return 0;
130 }
131 
132 /*TEST
133 
134    test:
135 
136    test:
137       suffix: 2
138       args: -bs 2
139 
140 TEST*/
141