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