xref: /petsc/src/mat/tests/ex26.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQBAIJ,&isseqbaij));
15   CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQSBAIJ,&isseqsbaij));
16   if (isseqbaij || isseqsbaij) {
17     CHKERRQ(MatGetBlockSize(A,&bs));
18   }
19   CHKERRQ(MatGetType(A,&type));
20   CHKERRQ(MatGetRowIJ(A,shift,symmetric,compressed,&nr,&ia,&ja,&done));
21   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"===========================================================\n"));
22   CHKERRQ(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     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT ":",i+shift));
25     for (j=ia[i];j<ia[i+1];j++) {
26       CHKERRQ(PetscPrintf(PETSC_COMM_SELF," %" PetscInt_FMT,ja[j-shift]));
27     }
28     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n"));
29   }
30   CHKERRQ(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   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
41   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
42   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
43   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
44   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
45   CHKERRQ(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   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
51   CHKERRQ(MatSetSizes(A,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
52   CHKERRQ(MatSetBlockSize(A,bs));
53   CHKERRQ(MatSetType(A,MATSEQAIJ));
54   CHKERRQ(MatSetUp(A));
55   CHKERRQ(MatCreate(PETSC_COMM_SELF,&B));
56   CHKERRQ(MatSetSizes(B,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
57   CHKERRQ(MatSetBlockSize(B,bs));
58   CHKERRQ(MatSetType(B,MATSEQBAIJ));
59   CHKERRQ(MatSetUp(B));
60   CHKERRQ(MatCreate(PETSC_COMM_SELF,&C));
61   CHKERRQ(MatSetSizes(C,m*n,m*n,PETSC_DECIDE,PETSC_DECIDE));
62   CHKERRQ(MatSetBlockSize(C,bs));
63   CHKERRQ(MatSetType(C,MATSEQSBAIJ));
64   CHKERRQ(MatSetUp(C));
65   CHKERRQ(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         CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
75         CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
76         CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
77       }
78       J = Ii + n;
79       if (J<m*n) {
80         CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
81         CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
82         CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
83       }
84       J = Ii - 1;
85       if (J>=0)  {
86         CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
87         CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
88         CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
89       }
90       J = Ii + 1;
91       if (J<m*n) {
92         CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
93         CHKERRQ(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES));
94         CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
95       }
96       v = 4.0;
97       CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
98       CHKERRQ(MatSetValues(B,1,&Ii,1,&Ii,&v,INSERT_VALUES));
99       CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
100     }
101   }
102   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
103   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
104   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
105   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
106   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
107   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
108 
109   /* test MatGetRowIJ for the three Mat types */
110   CHKERRQ(MatView(A,NULL));
111   CHKERRQ(MatView(B,NULL));
112   CHKERRQ(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         CHKERRQ(DumpCSR(A,shift,symmetric,compressed));
120         CHKERRQ(DumpCSR(B,shift,symmetric,compressed));
121         CHKERRQ(DumpCSR(C,shift,symmetric,compressed));
122       }
123     }
124   }
125   CHKERRQ(MatDestroy(&A));
126   CHKERRQ(MatDestroy(&B));
127   CHKERRQ(MatDestroy(&C));
128   CHKERRQ(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