xref: /petsc/src/mat/tests/ex99.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static const char help[] = "Tests MatCreateSubMatrix with MatSubMatrix versus MatAIJ, square, shifted (copied from ex97)\n";
2 
3 #include <petscmat.h>
4 
5 static PetscErrorCode AssembleMatrix(MPI_Comm comm, Mat *A) {
6   Mat      B;
7   PetscInt i, ms, me;
8 
9   PetscFunctionBegin;
10   PetscCall(MatCreate(comm, &B));
11   PetscCall(MatSetSizes(B, 6, 6, PETSC_DETERMINE, PETSC_DETERMINE));
12   PetscCall(MatSetFromOptions(B));
13   PetscCall(MatSetUp(B));
14   PetscCall(MatGetOwnershipRange(B, &ms, &me));
15   for (i = ms; i < me; i++) { PetscCall(MatSetValue(B, i, i, 1.0 * i, INSERT_VALUES)); }
16   PetscCall(MatSetValue(B, me - 1, me - 1, me * me, INSERT_VALUES));
17   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
18   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
19   *A = B;
20   PetscFunctionReturn(0);
21 }
22 
23 static PetscErrorCode Compare2(Vec *X, const char *test) {
24   PetscReal norm;
25   Vec       Y;
26   PetscInt  verbose = 0;
27 
28   PetscFunctionBegin;
29   PetscCall(VecDuplicate(X[0], &Y));
30   PetscCall(VecCopy(X[0], Y));
31   PetscCall(VecAYPX(Y, -1.0, X[1]));
32   PetscCall(VecNorm(Y, NORM_INFINITY, &norm));
33 
34   PetscCall(PetscOptionsGetInt(NULL, NULL, "-verbose", &verbose, NULL));
35   if (norm < PETSC_SQRT_MACHINE_EPSILON && verbose < 1) {
36     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%30s: norm difference < sqrt(eps_machine)\n", test));
37   } else {
38     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%30s: norm difference %g\n", test, (double)norm));
39   }
40   if (verbose > 1) {
41     PetscCall(VecView(X[0], PETSC_VIEWER_STDOUT_WORLD));
42     PetscCall(VecView(X[1], PETSC_VIEWER_STDOUT_WORLD));
43     PetscCall(VecView(Y, PETSC_VIEWER_STDOUT_WORLD));
44   }
45   PetscCall(VecDestroy(&Y));
46   PetscFunctionReturn(0);
47 }
48 
49 static PetscErrorCode CheckMatrices(Mat A, Mat B, Vec left, Vec right, Vec X, Vec Y, Vec X1, Vec Y1) {
50   Vec *ltmp, *rtmp;
51 
52   PetscFunctionBegin;
53   PetscCall(VecDuplicateVecs(right, 2, &rtmp));
54   PetscCall(VecDuplicateVecs(left, 2, &ltmp));
55   PetscCall(MatScale(A, PETSC_PI));
56   PetscCall(MatScale(B, PETSC_PI));
57   PetscCall(MatDiagonalScale(A, left, right));
58   PetscCall(MatDiagonalScale(B, left, right));
59   PetscCall(MatShift(A, PETSC_PI));
60   PetscCall(MatShift(B, PETSC_PI));
61 
62   PetscCall(MatMult(A, X, ltmp[0]));
63   PetscCall(MatMult(B, X, ltmp[1]));
64   PetscCall(Compare2(ltmp, "MatMult"));
65 
66   PetscCall(MatMultTranspose(A, Y, rtmp[0]));
67   PetscCall(MatMultTranspose(B, Y, rtmp[1]));
68   PetscCall(Compare2(rtmp, "MatMultTranspose"));
69 
70   PetscCall(VecCopy(Y1, ltmp[0]));
71   PetscCall(VecCopy(Y1, ltmp[1]));
72   PetscCall(MatMultAdd(A, X, ltmp[0], ltmp[0]));
73   PetscCall(MatMultAdd(B, X, ltmp[1], ltmp[1]));
74   PetscCall(Compare2(ltmp, "MatMultAdd v2==v3"));
75 
76   PetscCall(MatMultAdd(A, X, Y1, ltmp[0]));
77   PetscCall(MatMultAdd(B, X, Y1, ltmp[1]));
78   PetscCall(Compare2(ltmp, "MatMultAdd v2!=v3"));
79 
80   PetscCall(VecCopy(X1, rtmp[0]));
81   PetscCall(VecCopy(X1, rtmp[1]));
82   PetscCall(MatMultTransposeAdd(A, Y, rtmp[0], rtmp[0]));
83   PetscCall(MatMultTransposeAdd(B, Y, rtmp[1], rtmp[1]));
84   PetscCall(Compare2(rtmp, "MatMultTransposeAdd v2==v3"));
85 
86   PetscCall(MatMultTransposeAdd(A, Y, X1, rtmp[0]));
87   PetscCall(MatMultTransposeAdd(B, Y, X1, rtmp[1]));
88   PetscCall(Compare2(rtmp, "MatMultTransposeAdd v2!=v3"));
89 
90   PetscCall(VecDestroyVecs(2, &ltmp));
91   PetscCall(VecDestroyVecs(2, &rtmp));
92   PetscFunctionReturn(0);
93 }
94 
95 int main(int argc, char *argv[]) {
96   Mat       A, B, Asub, Bsub;
97   PetscInt  ms, idxrow[3], idxcol[3];
98   Vec       left, right, X, Y, X1, Y1;
99   IS        isrow, iscol;
100   PetscBool random = PETSC_TRUE;
101 
102   PetscFunctionBeginUser;
103   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
104   PetscCall(AssembleMatrix(PETSC_COMM_WORLD, &A));
105   PetscCall(AssembleMatrix(PETSC_COMM_WORLD, &B));
106   PetscCall(MatSetOperation(B, MATOP_CREATE_SUBMATRIX, NULL));
107   PetscCall(MatSetOperation(B, MATOP_CREATE_SUBMATRICES, NULL));
108   PetscCall(MatGetOwnershipRange(A, &ms, NULL));
109 
110   idxrow[0] = ms + 1;
111   idxrow[1] = ms + 2;
112   idxrow[2] = ms + 4;
113   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 3, idxrow, PETSC_USE_POINTER, &isrow));
114 
115   idxcol[0] = ms + 1;
116   idxcol[1] = ms + 2;
117   idxcol[2] = ms + 4;
118   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 3, idxcol, PETSC_USE_POINTER, &iscol));
119 
120   PetscCall(MatCreateSubMatrix(A, isrow, iscol, MAT_INITIAL_MATRIX, &Asub));
121   PetscCall(MatCreateSubMatrix(B, isrow, iscol, MAT_INITIAL_MATRIX, &Bsub));
122 
123   PetscCall(MatCreateVecs(Asub, &right, &left));
124   PetscCall(VecDuplicate(right, &X));
125   PetscCall(VecDuplicate(right, &X1));
126   PetscCall(VecDuplicate(left, &Y));
127   PetscCall(VecDuplicate(left, &Y1));
128 
129   PetscCall(PetscOptionsGetBool(NULL, NULL, "-random", &random, NULL));
130   if (random) {
131     PetscCall(VecSetRandom(right, NULL));
132     PetscCall(VecSetRandom(left, NULL));
133     PetscCall(VecSetRandom(X, NULL));
134     PetscCall(VecSetRandom(Y, NULL));
135     PetscCall(VecSetRandom(X1, NULL));
136     PetscCall(VecSetRandom(Y1, NULL));
137   } else {
138     PetscCall(VecSet(right, 1.0));
139     PetscCall(VecSet(left, 2.0));
140     PetscCall(VecSet(X, 3.0));
141     PetscCall(VecSet(Y, 4.0));
142     PetscCall(VecSet(X1, 3.0));
143     PetscCall(VecSet(Y1, 4.0));
144   }
145   PetscCall(CheckMatrices(Asub, Bsub, left, right, X, Y, X1, Y1));
146   PetscCall(ISDestroy(&isrow));
147   PetscCall(ISDestroy(&iscol));
148   PetscCall(MatDestroy(&A));
149   PetscCall(MatDestroy(&B));
150   PetscCall(MatDestroy(&Asub));
151   PetscCall(MatDestroy(&Bsub));
152   PetscCall(VecDestroy(&left));
153   PetscCall(VecDestroy(&right));
154   PetscCall(VecDestroy(&X));
155   PetscCall(VecDestroy(&Y));
156   PetscCall(VecDestroy(&X1));
157   PetscCall(VecDestroy(&Y1));
158   PetscCall(PetscFinalize());
159   return 0;
160 }
161 
162 /*TEST
163 
164    test:
165       nsize: 3
166 
167 TEST*/
168