1*345a4b08SToby Isaac const char help[] = "Test MATDIAGONAL";
2*345a4b08SToby Isaac
3*345a4b08SToby Isaac #include <petsc/private/petscimpl.h>
4*345a4b08SToby Isaac #include <petscmat.h>
5*345a4b08SToby Isaac
main(int argc,char ** argv)6*345a4b08SToby Isaac int main(int argc, char **argv)
7*345a4b08SToby Isaac {
8*345a4b08SToby Isaac Vec a, a2, b, b2, c, c2, A_diag, A_inv_diag;
9*345a4b08SToby Isaac Mat A, B;
10*345a4b08SToby Isaac PetscInt n = 10;
11*345a4b08SToby Isaac
12*345a4b08SToby Isaac PetscCall(PetscInitialize(&argc, &argv, NULL, help));
13*345a4b08SToby Isaac
14*345a4b08SToby Isaac MPI_Comm comm = PETSC_COMM_SELF;
15*345a4b08SToby Isaac PetscCall(VecCreateSeq(comm, n, &a));
16*345a4b08SToby Isaac PetscCall(VecDuplicate(a, &b));
17*345a4b08SToby Isaac PetscCall(VecDuplicate(a, &c));
18*345a4b08SToby Isaac PetscRandom rand;
19*345a4b08SToby Isaac
20*345a4b08SToby Isaac PetscCall(PetscRandomCreate(comm, &rand));
21*345a4b08SToby Isaac PetscCall(VecSetRandom(a, rand));
22*345a4b08SToby Isaac PetscCall(VecSetRandom(b, rand));
23*345a4b08SToby Isaac
24*345a4b08SToby Isaac PetscCall(VecDuplicate(a, &a2));
25*345a4b08SToby Isaac PetscCall(VecCopy(a, a2));
26*345a4b08SToby Isaac PetscCall(VecDuplicate(b, &b2));
27*345a4b08SToby Isaac PetscCall(VecCopy(b, b2));
28*345a4b08SToby Isaac PetscCall(VecDuplicate(c, &c2));
29*345a4b08SToby Isaac
30*345a4b08SToby Isaac PetscCall(MatCreateDiagonal(a2, &A));
31*345a4b08SToby Isaac PetscCall(MatCreateDiagonal(b2, &B));
32*345a4b08SToby Isaac PetscCall(VecDestroy(&a2));
33*345a4b08SToby Isaac PetscCall(VecDestroy(&b2));
34*345a4b08SToby Isaac
35*345a4b08SToby Isaac PetscCall(VecDuplicate(a, &a2));
36*345a4b08SToby Isaac PetscCall(VecDuplicate(b, &b2));
37*345a4b08SToby Isaac
38*345a4b08SToby Isaac PetscCall(MatAXPY(A, 0.5, B, SAME_NONZERO_PATTERN));
39*345a4b08SToby Isaac PetscCall(VecAXPY(a, 0.5, b));
40*345a4b08SToby Isaac
41*345a4b08SToby Isaac PetscReal mat_norm, vec_norm;
42*345a4b08SToby Isaac PetscCall(VecNorm(a, NORM_2, &vec_norm));
43*345a4b08SToby Isaac PetscCall(MatNorm(A, NORM_FROBENIUS, &mat_norm));
44*345a4b08SToby Isaac PetscCheck(vec_norm == mat_norm, comm, PETSC_ERR_PLIB, "Norms don't match");
45*345a4b08SToby Isaac
46*345a4b08SToby Isaac // For diagonal matrix, all operator norms are the max norm of the vector
47*345a4b08SToby Isaac PetscCall(VecNorm(a, NORM_INFINITY, &vec_norm));
48*345a4b08SToby Isaac PetscCall(MatNorm(A, NORM_INFINITY, &mat_norm));
49*345a4b08SToby Isaac PetscCheck(vec_norm == mat_norm, comm, PETSC_ERR_PLIB, "Norms don't match");
50*345a4b08SToby Isaac PetscCall(MatNorm(A, NORM_1, &mat_norm));
51*345a4b08SToby Isaac PetscCheck(vec_norm == mat_norm, comm, PETSC_ERR_PLIB, "Norms don't match");
52*345a4b08SToby Isaac
53*345a4b08SToby Isaac PetscCall(VecPointwiseMult(c, b, a));
54*345a4b08SToby Isaac PetscCall(MatMult(A, b, c2));
55*345a4b08SToby Isaac PetscCall(VecAXPY(c2, -1.0, c));
56*345a4b08SToby Isaac PetscCall(VecNorm(c2, NORM_INFINITY, &vec_norm));
57*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatMult is not like VecPointwiseMultiply");
58*345a4b08SToby Isaac
59*345a4b08SToby Isaac PetscCall(VecPointwiseMult(c, b, a));
60*345a4b08SToby Isaac PetscCall(VecAXPY(c, 1.0, a));
61*345a4b08SToby Isaac PetscCall(MatMultAdd(A, b, a, c2));
62*345a4b08SToby Isaac PetscCall(VecAXPY(c2, -1.0, c));
63*345a4b08SToby Isaac PetscCall(VecNorm(c2, NORM_INFINITY, &vec_norm));
64*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatMultAdd gave unexpected value");
65*345a4b08SToby Isaac
66*345a4b08SToby Isaac PetscCall(VecSet(c, 1.2));
67*345a4b08SToby Isaac PetscCall(VecSet(c2, 1.2));
68*345a4b08SToby Isaac PetscCall(VecPointwiseMult(c, b, a));
69*345a4b08SToby Isaac PetscCall(VecAXPY(c, 1.0, c2));
70*345a4b08SToby Isaac PetscCall(MatMultAdd(A, b, c2, c2));
71*345a4b08SToby Isaac PetscCall(VecAXPY(c2, -1.0, c));
72*345a4b08SToby Isaac PetscCall(VecNorm(c2, NORM_INFINITY, &vec_norm));
73*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatMultAdd gave unexpected value");
74*345a4b08SToby Isaac
75*345a4b08SToby Isaac PetscCall(VecPointwiseDivide(c, b, a));
76*345a4b08SToby Isaac PetscCall(MatSolve(A, b, c2));
77*345a4b08SToby Isaac PetscCall(VecAXPY(c2, -1.0, c));
78*345a4b08SToby Isaac PetscCall(VecNorm(c2, NORM_INFINITY, &vec_norm));
79*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatMult is not like VecPointwiseMultiply");
80*345a4b08SToby Isaac
81*345a4b08SToby Isaac Mat A_dup;
82*345a4b08SToby Isaac PetscCall(MatDuplicate(A, MAT_DO_NOT_COPY_VALUES, &A_dup));
83*345a4b08SToby Isaac PetscCall(MatDestroy(&A_dup));
84*345a4b08SToby Isaac PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &A_dup));
85*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A_dup, a2));
86*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
87*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
88*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatDuplicate with MAT_COPY_VALUES did not make a duplicate vector");
89*345a4b08SToby Isaac PetscCall(MatDestroy(&A_dup));
90*345a4b08SToby Isaac
91*345a4b08SToby Isaac PetscCall(MatShift(A, 1.5));
92*345a4b08SToby Isaac PetscCall(VecShift(a, 1.5));
93*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A, a2));
94*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
95*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
96*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatShift gave different result from VecShift");
97*345a4b08SToby Isaac
98*345a4b08SToby Isaac PetscCall(MatScale(A, 0.75));
99*345a4b08SToby Isaac PetscCall(VecScale(a, 0.75));
100*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A, a2));
101*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
102*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
103*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatScale gave different result from VecScale");
104*345a4b08SToby Isaac
105*345a4b08SToby Isaac PetscCall(VecPointwiseMult(a, a, b));
106*345a4b08SToby Isaac PetscCall(MatDiagonalScale(A, b, NULL));
107*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A, a2));
108*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
109*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
110*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatDiagonalScale gave unexpected result");
111*345a4b08SToby Isaac
112*345a4b08SToby Isaac PetscCall(VecPointwiseMult(a, a, b));
113*345a4b08SToby Isaac PetscCall(MatDiagonalScale(A, NULL, b));
114*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A, a2));
115*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
116*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
117*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatDiagonalScale gave unexpected result");
118*345a4b08SToby Isaac
119*345a4b08SToby Isaac PetscCall(VecCopy(b, a));
120*345a4b08SToby Isaac PetscCall(MatDiagonalSet(A, b, INSERT_VALUES));
121*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A, a2));
122*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
123*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
124*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatDiagonalSet gave unexpected result");
125*345a4b08SToby Isaac
126*345a4b08SToby Isaac PetscCall(VecSetRandom(a, rand));
127*345a4b08SToby Isaac PetscCall(VecSetRandom(b, rand));
128*345a4b08SToby Isaac PetscCall(MatDiagonalSet(A, a, INSERT_VALUES));
129*345a4b08SToby Isaac PetscCall(VecAXPY(a, 1.0, b));
130*345a4b08SToby Isaac PetscCall(MatDiagonalSet(A, b, ADD_VALUES));
131*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A, a2));
132*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
133*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
134*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatDiagonalSet gave unexpected result");
135*345a4b08SToby Isaac
136*345a4b08SToby Isaac PetscCall(VecSetRandom(a, rand));
137*345a4b08SToby Isaac PetscCall(VecSetRandom(b, rand));
138*345a4b08SToby Isaac PetscCall(MatDiagonalSet(A, a, INSERT_VALUES));
139*345a4b08SToby Isaac PetscCall(VecPointwiseMax(a, a, b));
140*345a4b08SToby Isaac PetscCall(MatDiagonalSet(A, b, MAX_VALUES));
141*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A, a2));
142*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
143*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
144*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatDiagonalSet gave unexpected result");
145*345a4b08SToby Isaac
146*345a4b08SToby Isaac PetscCall(VecSetRandom(a, rand));
147*345a4b08SToby Isaac PetscCall(VecSetRandom(b, rand));
148*345a4b08SToby Isaac PetscCall(MatDiagonalSet(A, a, INSERT_VALUES));
149*345a4b08SToby Isaac PetscCall(VecPointwiseMin(a, a, b));
150*345a4b08SToby Isaac PetscCall(MatDiagonalSet(A, b, MIN_VALUES));
151*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A, a2));
152*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
153*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
154*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatDiagonalSet gave unexpected result");
155*345a4b08SToby Isaac
156*345a4b08SToby Isaac PetscCall(VecSetRandom(a, rand));
157*345a4b08SToby Isaac PetscCall(VecSetRandom(b, rand));
158*345a4b08SToby Isaac PetscCall(MatDiagonalSet(A, a, INSERT_VALUES));
159*345a4b08SToby Isaac PetscCall(MatDiagonalSet(A, b, NOT_SET_VALUES));
160*345a4b08SToby Isaac PetscCall(MatGetDiagonal(A, a2));
161*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, a));
162*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
163*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatDiagonalSet gave unexpected result");
164*345a4b08SToby Isaac
165*345a4b08SToby Isaac PetscCall(VecSet(a2, 0.5));
166*345a4b08SToby Isaac
167*345a4b08SToby Isaac PetscObjectState state_pre, state_post;
168*345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)A, &state_pre));
169*345a4b08SToby Isaac PetscCall(MatDiagonalGetInverseDiagonal(A, &A_inv_diag));
170*345a4b08SToby Isaac PetscCall(MatDiagonalRestoreInverseDiagonal(A, &A_inv_diag));
171*345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)A, &state_post));
172*345a4b08SToby Isaac PetscCheck(state_pre == state_post, comm, PETSC_ERR_PLIB, "State changed on noop");
173*345a4b08SToby Isaac
174*345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)A, &state_pre));
175*345a4b08SToby Isaac PetscCall(MatDiagonalGetInverseDiagonal(A, &A_inv_diag));
176*345a4b08SToby Isaac PetscCall(VecSet(A_inv_diag, 2.0));
177*345a4b08SToby Isaac PetscCall(MatDiagonalRestoreInverseDiagonal(A, &A_inv_diag));
178*345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)A, &state_post));
179*345a4b08SToby Isaac PetscCheck(state_pre != state_post, comm, PETSC_ERR_PLIB, "State not changed on mutation");
180*345a4b08SToby Isaac
181*345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)A, &state_pre));
182*345a4b08SToby Isaac PetscCall(MatDiagonalGetDiagonal(A, &A_diag));
183*345a4b08SToby Isaac PetscCall(MatDiagonalRestoreDiagonal(A, &A_diag));
184*345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)A, &state_post));
185*345a4b08SToby Isaac PetscCheck(state_pre == state_post, comm, PETSC_ERR_PLIB, "State changed on noop");
186*345a4b08SToby Isaac
187*345a4b08SToby Isaac PetscCall(MatDiagonalGetDiagonal(A, &A_diag));
188*345a4b08SToby Isaac PetscCall(VecAXPY(a2, -1.0, A_diag));
189*345a4b08SToby Isaac PetscCall(VecSet(A_diag, 1.0));
190*345a4b08SToby Isaac PetscCall(MatDiagonalRestoreDiagonal(A, &A_diag));
191*345a4b08SToby Isaac PetscCall(PetscObjectStateGet((PetscObject)A, &state_post));
192*345a4b08SToby Isaac PetscCheck(state_pre != state_post, comm, PETSC_ERR_PLIB, "State not changed on mutation");
193*345a4b08SToby Isaac
194*345a4b08SToby Isaac PetscCall(VecNorm(a2, NORM_INFINITY, &vec_norm));
195*345a4b08SToby Isaac PetscCheck(vec_norm < PETSC_SMALL, comm, PETSC_ERR_PLIB, "MatDiagonalGetInverse gave unexpected result");
196*345a4b08SToby Isaac
197*345a4b08SToby Isaac PetscCall(MatZeroEntries(A));
198*345a4b08SToby Isaac PetscCall(MatNorm(A, NORM_INFINITY, &mat_norm));
199*345a4b08SToby Isaac PetscCheck(mat_norm == 0.0, comm, PETSC_ERR_PLIB, "MatZeroEntries gave unexpected result");
200*345a4b08SToby Isaac PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF));
201*345a4b08SToby Isaac PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_ASCII_INFO));
202*345a4b08SToby Isaac PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF));
203*345a4b08SToby Isaac PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF));
204*345a4b08SToby Isaac
205*345a4b08SToby Isaac PetscCall(MatDestroy(&A));
206*345a4b08SToby Isaac PetscCall(MatDestroy(&B));
207*345a4b08SToby Isaac
208*345a4b08SToby Isaac PetscCall(PetscRandomDestroy(&rand));
209*345a4b08SToby Isaac PetscCall(VecDestroy(&c2));
210*345a4b08SToby Isaac PetscCall(VecDestroy(&b2));
211*345a4b08SToby Isaac PetscCall(VecDestroy(&a2));
212*345a4b08SToby Isaac PetscCall(VecDestroy(&c));
213*345a4b08SToby Isaac PetscCall(VecDestroy(&b));
214*345a4b08SToby Isaac PetscCall(VecDestroy(&a));
215*345a4b08SToby Isaac PetscCall(PetscFinalize());
216*345a4b08SToby Isaac return 0;
217*345a4b08SToby Isaac }
218*345a4b08SToby Isaac
219*345a4b08SToby Isaac /*TEST
220*345a4b08SToby Isaac
221*345a4b08SToby Isaac test:
222*345a4b08SToby Isaac suffix: 0
223*345a4b08SToby Isaac
224*345a4b08SToby Isaac TEST*/
225