xref: /petsc/src/mat/tests/ex74.c (revision 9a3a8673b4aea812b2f0c314666d2e7ff14d2577)
1c4762a1bSJed Brown static char help[] = "Tests the various sequential routines in MATSEQSBAIJ format.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   PetscMPIInt   size;
8c4762a1bSJed Brown   Vec           x, y, b, s1, s2;
9c4762a1bSJed Brown   Mat           A;                     /* linear system matrix */
10c4762a1bSJed Brown   Mat           sA, sB, sFactor, B, C; /* symmetric matrices */
11c4762a1bSJed Brown   PetscInt      n, mbs = 16, bs = 1, nz = 3, prob = 1, i, j, k1, k2, col[3], lf, block, row, Ii, J, n1, inc;
12c4762a1bSJed Brown   PetscReal     norm1, norm2, rnorm, tol = 10 * PETSC_SMALL;
13c4762a1bSJed Brown   PetscScalar   neg_one = -1.0, four = 4.0, value[3];
14c4762a1bSJed Brown   IS            perm, iscol;
15c4762a1bSJed Brown   PetscRandom   rdm;
16c4762a1bSJed Brown   PetscBool     doIcc = PETSC_TRUE, equal;
17c4762a1bSJed Brown   MatInfo       minfo1, minfo2;
18c4762a1bSJed Brown   MatFactorInfo factinfo;
19c4762a1bSJed Brown   MatType       type;
20c4762a1bSJed Brown 
21327415f7SBarry Smith   PetscFunctionBeginUser;
22c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mbs", &mbs, NULL));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   n = mbs * bs;
299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
309566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
319566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATSEQBAIJ));
329566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
339566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(A, bs, nz, NULL));
34c4762a1bSJed Brown 
359566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &sA));
369566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(sA, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
379566063dSJacob Faibussowitsch   PetscCall(MatSetType(sA, MATSEQSBAIJ));
389566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(sA));
399566063dSJacob Faibussowitsch   PetscCall(MatGetType(sA, &type));
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sA, MATSEQSBAIJ, &doIcc));
419566063dSJacob Faibussowitsch   PetscCall(MatSeqSBAIJSetPreallocation(sA, bs, nz, NULL));
429566063dSJacob Faibussowitsch   PetscCall(MatSetOption(sA, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Test MatGetOwnershipRange() */
459566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &Ii, &J));
469566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(sA, &i, &j));
4748a46eb9SPierre Jolivet   if (i - Ii || j - J) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetOwnershipRange() in MatSBAIJ format\n"));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /* Assemble matrix */
50c4762a1bSJed Brown   if (bs == 1) {
519566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_problem", &prob, NULL));
52c4762a1bSJed Brown     if (prob == 1) { /* tridiagonal matrix */
539371c9d4SSatish Balay       value[0] = -1.0;
549371c9d4SSatish Balay       value[1] = 2.0;
559371c9d4SSatish Balay       value[2] = -1.0;
56c4762a1bSJed Brown       for (i = 1; i < n - 1; i++) {
579371c9d4SSatish Balay         col[0] = i - 1;
589371c9d4SSatish Balay         col[1] = i;
599371c9d4SSatish Balay         col[2] = i + 1;
609566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
619566063dSJacob Faibussowitsch         PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES));
62c4762a1bSJed Brown       }
639371c9d4SSatish Balay       i      = n - 1;
649371c9d4SSatish Balay       col[0] = 0;
659371c9d4SSatish Balay       col[1] = n - 2;
669371c9d4SSatish Balay       col[2] = n - 1;
67c4762a1bSJed Brown 
689371c9d4SSatish Balay       value[0] = 0.1;
699371c9d4SSatish Balay       value[1] = -1;
709371c9d4SSatish Balay       value[2] = 2;
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
739566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown       i        = 0;
769371c9d4SSatish Balay       col[0]   = n - 1;
779371c9d4SSatish Balay       col[1]   = 1;
789371c9d4SSatish Balay       col[2]   = 0;
799371c9d4SSatish Balay       value[0] = 0.1;
809371c9d4SSatish Balay       value[1] = -1.0;
819371c9d4SSatish Balay       value[2] = 2;
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
849566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown     } else if (prob == 2) { /* matrix for the five point stencil */
87c4762a1bSJed Brown       n1 = (PetscInt)(PetscSqrtReal((PetscReal)n) + 0.001);
88bc30f867SBarry Smith       PetscCheck(n1 * n1 == n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "sqrt(n) must be a positive integer!");
89c4762a1bSJed Brown       for (i = 0; i < n1; i++) {
90c4762a1bSJed Brown         for (j = 0; j < n1; j++) {
91c4762a1bSJed Brown           Ii = j + n1 * i;
92c4762a1bSJed Brown           if (i > 0) {
93c4762a1bSJed Brown             J = Ii - n1;
949566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
959566063dSJacob Faibussowitsch             PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
96c4762a1bSJed Brown           }
97c4762a1bSJed Brown           if (i < n1 - 1) {
98c4762a1bSJed Brown             J = Ii + n1;
999566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
1009566063dSJacob Faibussowitsch             PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
101c4762a1bSJed Brown           }
102c4762a1bSJed Brown           if (j > 0) {
103c4762a1bSJed Brown             J = Ii - 1;
1049566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
1059566063dSJacob Faibussowitsch             PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
106c4762a1bSJed Brown           }
107c4762a1bSJed Brown           if (j < n1 - 1) {
108c4762a1bSJed Brown             J = Ii + 1;
1099566063dSJacob Faibussowitsch             PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
1109566063dSJacob Faibussowitsch             PetscCall(MatSetValues(sA, 1, &Ii, 1, &J, &neg_one, INSERT_VALUES));
111c4762a1bSJed Brown           }
1129566063dSJacob Faibussowitsch           PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &four, INSERT_VALUES));
1139566063dSJacob Faibussowitsch           PetscCall(MatSetValues(sA, 1, &Ii, 1, &Ii, &four, INSERT_VALUES));
114c4762a1bSJed Brown         }
115c4762a1bSJed Brown       }
116c4762a1bSJed Brown     }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   } else { /* bs > 1 */
119c4762a1bSJed Brown     for (block = 0; block < n / bs; block++) {
120c4762a1bSJed Brown       /* diagonal blocks */
1219371c9d4SSatish Balay       value[0] = -1.0;
1229371c9d4SSatish Balay       value[1] = 4.0;
1239371c9d4SSatish Balay       value[2] = -1.0;
124c4762a1bSJed Brown       for (i = 1 + block * bs; i < bs - 1 + block * bs; i++) {
1259371c9d4SSatish Balay         col[0] = i - 1;
1269371c9d4SSatish Balay         col[1] = i;
1279371c9d4SSatish Balay         col[2] = i + 1;
1289566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
1299566063dSJacob Faibussowitsch         PetscCall(MatSetValues(sA, 1, &i, 3, col, value, INSERT_VALUES));
130c4762a1bSJed Brown       }
1319371c9d4SSatish Balay       i      = bs - 1 + block * bs;
1329371c9d4SSatish Balay       col[0] = bs - 2 + block * bs;
1339371c9d4SSatish Balay       col[1] = bs - 1 + block * bs;
134c4762a1bSJed Brown 
1359371c9d4SSatish Balay       value[0] = -1.0;
1369371c9d4SSatish Balay       value[1] = 4.0;
137c4762a1bSJed Brown 
1389566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
1399566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA, 1, &i, 2, col, value, INSERT_VALUES));
140c4762a1bSJed Brown 
1419371c9d4SSatish Balay       i      = 0 + block * bs;
1429371c9d4SSatish Balay       col[0] = 0 + block * bs;
1439371c9d4SSatish Balay       col[1] = 1 + block * bs;
144c4762a1bSJed Brown 
1459371c9d4SSatish Balay       value[0] = 4.0;
1469371c9d4SSatish Balay       value[1] = -1.0;
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
1499566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA, 1, &i, 2, col, value, INSERT_VALUES));
150c4762a1bSJed Brown     }
151c4762a1bSJed Brown     /* off-diagonal blocks */
152c4762a1bSJed Brown     value[0] = -1.0;
153c4762a1bSJed Brown     for (i = 0; i < (n / bs - 1) * bs; i++) {
154c4762a1bSJed Brown       col[0] = i + bs;
155c4762a1bSJed Brown 
1569566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &i, 1, col, value, INSERT_VALUES));
1579566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA, 1, &i, 1, col, value, INSERT_VALUES));
158c4762a1bSJed Brown 
1599371c9d4SSatish Balay       col[0] = i;
1609371c9d4SSatish Balay       row    = i + bs;
161c4762a1bSJed Brown 
1629566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &row, 1, col, value, INSERT_VALUES));
1639566063dSJacob Faibussowitsch       PetscCall(MatSetValues(sA, 1, &row, 1, col, value, INSERT_VALUES));
164c4762a1bSJed Brown     }
165c4762a1bSJed Brown   }
1669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
168c4762a1bSJed Brown 
1699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(sA, MAT_FINAL_ASSEMBLY));
1709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(sA, MAT_FINAL_ASSEMBLY));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* Test MatGetInfo() of A and sA */
1739566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_LOCAL, &minfo1));
1749566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(sA, MAT_LOCAL, &minfo2));
175c4762a1bSJed Brown   i  = (int)(minfo1.nz_used - minfo2.nz_used);
176c4762a1bSJed Brown   j  = (int)(minfo1.nz_allocated - minfo2.nz_allocated);
177c4762a1bSJed Brown   k1 = (int)(minfo1.nz_allocated - minfo1.nz_used);
178c4762a1bSJed Brown   k2 = (int)(minfo2.nz_allocated - minfo2.nz_used);
17948a46eb9SPierre Jolivet   if (i < 0 || j < 0 || k1 < 0 || k2 < 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error (compare A and sA): MatGetInfo()\n"));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* Test MatDuplicate() */
1829566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_FROBENIUS, &norm1));
1839566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(sA, MAT_COPY_VALUES, &sB));
1849566063dSJacob Faibussowitsch   PetscCall(MatEqual(sA, sB, &equal));
18528b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in MatDuplicate()");
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* Test MatNorm() */
1889566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_FROBENIUS, &norm1));
1899566063dSJacob Faibussowitsch   PetscCall(MatNorm(sB, NORM_FROBENIUS, &norm2));
190c4762a1bSJed Brown   rnorm = PetscAbsReal(norm1 - norm2) / norm2;
19148a46eb9SPierre Jolivet   if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_FROBENIUS, NormA=%16.14e NormsB=%16.14e\n", (double)norm1, (double)norm2));
1929566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_INFINITY, &norm1));
1939566063dSJacob Faibussowitsch   PetscCall(MatNorm(sB, NORM_INFINITY, &norm2));
194c4762a1bSJed Brown   rnorm = PetscAbsReal(norm1 - norm2) / norm2;
19548a46eb9SPierre Jolivet   if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n", (double)norm1, (double)norm2));
1969566063dSJacob Faibussowitsch   PetscCall(MatNorm(A, NORM_1, &norm1));
1979566063dSJacob Faibussowitsch   PetscCall(MatNorm(sB, NORM_1, &norm2));
198c4762a1bSJed Brown   rnorm = PetscAbsReal(norm1 - norm2) / norm2;
19948a46eb9SPierre Jolivet   if (rnorm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatNorm_INFINITY(), NormA=%16.14e NormsB=%16.14e\n", (double)norm1, (double)norm2));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */
2029566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(A, MAT_LOCAL, &minfo1));
2039566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(sB, MAT_LOCAL, &minfo2));
204c4762a1bSJed Brown   i  = (int)(minfo1.nz_used - minfo2.nz_used);
205c4762a1bSJed Brown   j  = (int)(minfo1.nz_allocated - minfo2.nz_allocated);
206c4762a1bSJed Brown   k1 = (int)(minfo1.nz_allocated - minfo1.nz_used);
207c4762a1bSJed Brown   k2 = (int)(minfo2.nz_allocated - minfo2.nz_used);
20848a46eb9SPierre Jolivet   if (i < 0 || j < 0 || k1 < 0 || k2 < 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error(compare A and sB): MatGetInfo()\n"));
209c4762a1bSJed Brown 
2109566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &Ii, &J));
2119566063dSJacob Faibussowitsch   PetscCall(MatGetSize(sB, &i, &j));
21248a46eb9SPierre Jolivet   if (i - Ii || j - J) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetSize()\n"));
213c4762a1bSJed Brown 
2149566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &Ii));
2159566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(sB, &i));
21648a46eb9SPierre Jolivet   if (i - Ii) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatGetBlockSize()\n"));
217c4762a1bSJed Brown 
2189566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
2199566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rdm));
2209566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &x));
2219566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &s1));
2229566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &s2));
2239566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &y));
2249566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &b));
2259566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(x, rdm));
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */
228c4762a1bSJed Brown #if !defined(PETSC_USE_COMPLEX)
229c4762a1bSJed Brown   /* Scaling matrix with complex numbers results non-spd matrix,
230c4762a1bSJed Brown      causing crash of MatForwardSolve() and MatBackwardSolve() */
2319566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(A, x, x));
2329566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(sB, x, x));
2339566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A, sB, 10, &equal));
23428b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Error in MatDiagonalScale");
235c4762a1bSJed Brown 
2369566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(A, s1));
2379566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(sB, s2));
2389566063dSJacob Faibussowitsch   PetscCall(VecAXPY(s2, neg_one, s1));
2399566063dSJacob Faibussowitsch   PetscCall(VecNorm(s2, NORM_1, &norm1));
24048a46eb9SPierre Jolivet   if (norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatGetDiagonal(), ||s1-s2||=%g\n", (double)norm1));
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   {
243c4762a1bSJed Brown     PetscScalar alpha = 0.1;
2449566063dSJacob Faibussowitsch     PetscCall(MatScale(A, alpha));
2459566063dSJacob Faibussowitsch     PetscCall(MatScale(sB, alpha));
246c4762a1bSJed Brown   }
247c4762a1bSJed Brown #endif
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* Test MatGetRowMaxAbs() */
2509566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(A, s1, NULL));
2519566063dSJacob Faibussowitsch   PetscCall(MatGetRowMaxAbs(sB, s2, NULL));
2529566063dSJacob Faibussowitsch   PetscCall(VecNorm(s1, NORM_1, &norm1));
2539566063dSJacob Faibussowitsch   PetscCall(VecNorm(s2, NORM_1, &norm2));
254c4762a1bSJed Brown   norm1 -= norm2;
25548a46eb9SPierre Jolivet   if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatGetRowMaxAbs() \n"));
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /* Test MatMult() */
258c4762a1bSJed Brown   for (i = 0; i < 40; i++) {
2599566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rdm));
2609566063dSJacob Faibussowitsch     PetscCall(MatMult(A, x, s1));
2619566063dSJacob Faibussowitsch     PetscCall(MatMult(sB, x, s2));
2629566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_1, &norm1));
2639566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_1, &norm2));
264c4762a1bSJed Brown     norm1 -= norm2;
26548a46eb9SPierre Jolivet     if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error: MatMult(), norm1-norm2: %g\n", (double)norm1));
266c4762a1bSJed Brown   }
267c4762a1bSJed Brown 
268c4762a1bSJed Brown   /* MatMultAdd() */
269c4762a1bSJed Brown   for (i = 0; i < 40; i++) {
2709566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(x, rdm));
2719566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(y, rdm));
2729566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(A, x, y, s1));
2739566063dSJacob Faibussowitsch     PetscCall(MatMultAdd(sB, x, y, s2));
2749566063dSJacob Faibussowitsch     PetscCall(VecNorm(s1, NORM_1, &norm1));
2759566063dSJacob Faibussowitsch     PetscCall(VecNorm(s2, NORM_1, &norm2));
276c4762a1bSJed Brown     norm1 -= norm2;
27748a46eb9SPierre Jolivet     if (norm1 < -tol || norm1 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatMultAdd(), norm1-norm2: %g\n", (double)norm1));
278c4762a1bSJed Brown   }
279c4762a1bSJed Brown 
280c20d7725SJed Brown   /* Test MatMatMult() for sbaij and dense matrices */
2819566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, n, 5 * n, NULL, &B));
2829566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(B, rdm));
283fb842aefSJose E. Roman   PetscCall(MatMatMult(sA, B, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &C));
2849566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(sA, B, C, 5 * n, &equal));
28528b400f6SJacob Faibussowitsch   PetscCheck(equal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error: MatMatMult()");
2869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
2879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */
2909566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A, MATORDERINGNATURAL, &perm, &iscol));
2919566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
292c4762a1bSJed Brown   norm1 = tol;
293c4762a1bSJed Brown   inc   = bs;
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   /* initialize factinfo */
2969566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&factinfo, sizeof(MatFactorInfo)));
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   for (lf = -1; lf < 10; lf += inc) {
299c4762a1bSJed Brown     if (lf == -1) { /* Cholesky factor of sB (duplicate sA) */
300c4762a1bSJed Brown       factinfo.fill = 5.0;
301c4762a1bSJed Brown 
3029566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sB, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &sFactor));
3039566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(sFactor, sB, perm, &factinfo));
304c4762a1bSJed Brown     } else if (!doIcc) break;
3059371c9d4SSatish Balay     else { /* incomplete Cholesky factor */ factinfo.fill = 5.0;
306c4762a1bSJed Brown       factinfo.levels                                     = lf;
307c4762a1bSJed Brown 
3089566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(sB, MATSOLVERPETSC, MAT_FACTOR_ICC, &sFactor));
3099566063dSJacob Faibussowitsch       PetscCall(MatICCFactorSymbolic(sFactor, sB, perm, &factinfo));
310c4762a1bSJed Brown     }
3119566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(sFactor, sB, &factinfo));
312c4762a1bSJed Brown     /* MatView(sFactor, PETSC_VIEWER_DRAW_WORLD); */
313c4762a1bSJed Brown 
314c4762a1bSJed Brown     /* test MatGetDiagonal on numeric factor */
315c4762a1bSJed Brown     /*
316c4762a1bSJed Brown     if (lf == -1) {
3179566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(sFactor,s1));
318c4762a1bSJed Brown       printf(" in ex74.c, diag: \n");
3199566063dSJacob Faibussowitsch       PetscCall(VecView(s1,PETSC_VIEWER_STDOUT_SELF));
320c4762a1bSJed Brown     }
321c4762a1bSJed Brown     */
322c4762a1bSJed Brown 
3239566063dSJacob Faibussowitsch     PetscCall(MatMult(sB, x, b));
324c4762a1bSJed Brown 
325c4762a1bSJed Brown     /* test MatForwardSolve() and MatBackwardSolve() */
326c4762a1bSJed Brown     if (lf == -1) {
3279566063dSJacob Faibussowitsch       PetscCall(MatForwardSolve(sFactor, b, s1));
3289566063dSJacob Faibussowitsch       PetscCall(MatBackwardSolve(sFactor, s1, s2));
3299566063dSJacob Faibussowitsch       PetscCall(VecAXPY(s2, neg_one, x));
3309566063dSJacob Faibussowitsch       PetscCall(VecNorm(s2, NORM_2, &norm2));
33148a46eb9SPierre Jolivet       if (10 * norm1 < norm2) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatForwardSolve and BackwardSolve: Norm of error=%g, bs=%" PetscInt_FMT "\n", (double)norm2, bs));
332c4762a1bSJed Brown     }
333c4762a1bSJed Brown 
334c4762a1bSJed Brown     /* test MatSolve() */
3359566063dSJacob Faibussowitsch     PetscCall(MatSolve(sFactor, b, y));
3369566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&sFactor));
337c4762a1bSJed Brown     /* Check the error */
3389566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, neg_one, x));
3399566063dSJacob Faibussowitsch     PetscCall(VecNorm(y, NORM_2, &norm2));
34048a46eb9SPierre Jolivet     if (10 * norm1 < norm2 && lf - inc != -1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "lf=%" PetscInt_FMT ", %" PetscInt_FMT ", Norm of error=%g, %g\n", lf - inc, lf, (double)norm1, (double)norm2));
341c4762a1bSJed Brown     norm1 = norm2;
342c4762a1bSJed Brown     if (norm2 < tol && lf != -1) break;
343c4762a1bSJed Brown   }
344c4762a1bSJed Brown 
345c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
346*c13ae4d5SJunchao Zhang 
347*c13ae4d5SJunchao Zhang   #if defined(PETSC_USE_REAL___FLOAT128)
348*c13ae4d5SJunchao Zhang   tol = 1e-10; // since MUMPS is run in double
349*c13ae4d5SJunchao Zhang   #endif
350*c13ae4d5SJunchao Zhang 
3519566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(sA, MATSOLVERMUMPS, MAT_FACTOR_CHOLESKY, &sFactor));
3529566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorSymbolic(sFactor, sA, NULL, NULL));
3539566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorNumeric(sFactor, sA, NULL));
354c4762a1bSJed Brown   for (i = 0; i < 10; i++) {
3559566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(b, rdm));
3569566063dSJacob Faibussowitsch     PetscCall(MatSolve(sFactor, b, y));
357c4762a1bSJed Brown     /* Check the error */
3589566063dSJacob Faibussowitsch     PetscCall(MatMult(sA, y, x));
3599566063dSJacob Faibussowitsch     PetscCall(VecAXPY(x, neg_one, b));
3609566063dSJacob Faibussowitsch     PetscCall(VecNorm(x, NORM_2, &norm2));
36148a46eb9SPierre Jolivet     if (norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Error:MatSolve(), norm2: %g\n", (double)norm2));
362c4762a1bSJed Brown   }
3639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sFactor));
364c4762a1bSJed Brown #endif
365c4762a1bSJed Brown 
3669566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
367c4762a1bSJed Brown 
3689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sB));
3709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
3719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
3739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s1));
3749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&s2));
3759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
3769566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rdm));
377c4762a1bSJed Brown 
3789566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
379b122ec5aSJacob Faibussowitsch   return 0;
380c4762a1bSJed Brown }
381c4762a1bSJed Brown 
382c4762a1bSJed Brown /*TEST
383c4762a1bSJed Brown 
384c4762a1bSJed Brown    test:
385c4762a1bSJed Brown       args: -bs {{1 2 3 4 5 6 7 8}}
3863886731fSPierre Jolivet       output_file: output/empty.out
387c4762a1bSJed Brown 
388c4762a1bSJed Brown TEST*/
389