xref: /petsc/src/mat/tests/ex192.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Tests MatSolve() and MatMatSolve() with MUMPS or MKL_PARDISO sequential solvers in Schur complement mode.\n\
2c4762a1bSJed Brown Example: mpiexec -n 1 ./ex192 -f <matrix binary file> -nrhs 4 -symmetric_solve -hermitian_solve -schur_ratio 0.3\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
69371c9d4SSatish Balay int main(int argc, char **args) {
7c4762a1bSJed Brown   Mat         A, RHS, C, F, X, S;
8c4762a1bSJed Brown   Vec         u, x, b;
9c4762a1bSJed Brown   Vec         xschur, bschur, uschur;
10c4762a1bSJed Brown   IS          is_schur;
11c4762a1bSJed Brown   PetscMPIInt size;
12c4762a1bSJed Brown   PetscInt    isolver = 0, size_schur, m, n, nfact, nsolve, nrhs;
13c4762a1bSJed Brown   PetscReal   norm, tol = PETSC_SQRT_MACHINE_EPSILON;
14c4762a1bSJed Brown   PetscRandom rand;
15c4762a1bSJed Brown   PetscBool   data_provided, herm, symm, use_lu, cuda = PETSC_FALSE;
16c4762a1bSJed Brown   PetscReal   sratio = 5.1 / 12.;
17c4762a1bSJed Brown   PetscViewer fd; /* viewer */
18c4762a1bSJed Brown   char        solver[256];
19c4762a1bSJed Brown   char        file[PETSC_MAX_PATH_LEN]; /* input file name */
20c4762a1bSJed Brown 
21327415f7SBarry Smith   PetscFunctionBeginUser;
229566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, 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 test");
25c4762a1bSJed Brown   /* Determine which type of solver we want to test for */
26c4762a1bSJed Brown   herm = PETSC_FALSE;
27c4762a1bSJed Brown   symm = PETSC_FALSE;
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-symmetric_solve", &symm, NULL));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-hermitian_solve", &herm, NULL));
30c4762a1bSJed Brown   if (herm) symm = PETSC_TRUE;
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-cuda_solve", &cuda, NULL));
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &data_provided));
36c4762a1bSJed Brown   if (!data_provided) { /* get matrices from PETSc distribution */
379566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(file, "${PETSC_DIR}/share/petsc/datafiles/matrices/", sizeof(file)));
38c4762a1bSJed Brown     if (symm) {
39c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
409566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(file, "hpd-complex-", sizeof(file)));
41c4762a1bSJed Brown #else
429566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(file, "spd-real-", sizeof(file)));
43c4762a1bSJed Brown #endif
44c4762a1bSJed Brown     } else {
45c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
469566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(file, "nh-complex-", sizeof(file)));
47c4762a1bSJed Brown #else
489566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(file, "ns-real-", sizeof(file)));
49c4762a1bSJed Brown #endif
50c4762a1bSJed Brown     }
51c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
529566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(file, "int64-", sizeof(file)));
53c4762a1bSJed Brown #else
549566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(file, "int32-", sizeof(file)));
55c4762a1bSJed Brown #endif
56c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
579566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(file, "float32", sizeof(file)));
58c4762a1bSJed Brown #else
599566063dSJacob Faibussowitsch     PetscCall(PetscStrlcat(file, "float64", sizeof(file)));
60c4762a1bSJed Brown #endif
61c4762a1bSJed Brown   }
62c4762a1bSJed Brown   /* Load matrix A */
639566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
649566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
659566063dSJacob Faibussowitsch   PetscCall(MatLoad(A, fd));
669566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
679566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &m, &n));
6808401ef6SPierre Jolivet   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n);
69c4762a1bSJed Brown 
70a5b23f4aSJose E. Roman   /* Create dense matrix C and X; C holds true solution with identical columns */
71c4762a1bSJed Brown   nrhs = 2;
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL));
739566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
749566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, m, PETSC_DECIDE, PETSC_DECIDE, nrhs));
759566063dSJacob Faibussowitsch   PetscCall(MatSetType(C, MATDENSE));
769566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
779566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
809566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
819566063dSJacob Faibussowitsch   PetscCall(MatSetRandom(C, rand));
829566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &X));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* Create vectors */
859566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
869566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
879566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
889566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &b));
899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &u)); /* save the true solution */
90c4762a1bSJed Brown 
919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-solver", &isolver, NULL));
92c4762a1bSJed Brown   switch (isolver) {
93c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
949371c9d4SSatish Balay   case 0: PetscCall(PetscStrcpy(solver, MATSOLVERMUMPS)); break;
95c4762a1bSJed Brown #endif
96c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_PARDISO)
979371c9d4SSatish Balay   case 1: PetscCall(PetscStrcpy(solver, MATSOLVERMKL_PARDISO)); break;
98c4762a1bSJed Brown #endif
999371c9d4SSatish Balay   default: PetscCall(PetscStrcpy(solver, MATSOLVERPETSC)); break;
100c4762a1bSJed Brown   }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
103c4762a1bSJed Brown   if (isolver == 0 && symm && !data_provided) { /* MUMPS (5.0.0) does not have support for hermitian matrices, so make them symmetric */
104c4762a1bSJed Brown     PetscScalar im  = PetscSqrtScalar((PetscScalar)-1.);
105c4762a1bSJed Brown     PetscScalar val = -1.0;
106c4762a1bSJed Brown     val             = val + im;
1079566063dSJacob Faibussowitsch     PetscCall(MatSetValue(A, 1, 0, val, INSERT_VALUES));
1089566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1099566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown #endif
112c4762a1bSJed Brown 
1139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-schur_ratio", &sratio, NULL));
114bc30f867SBarry Smith   PetscCheck(sratio >= 0. && sratio <= 1., PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid ratio for schur degrees of freedom %g", (double)sratio);
115c4762a1bSJed Brown   size_schur = (PetscInt)(sratio * m);
116c4762a1bSJed Brown 
1179566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Solving with %s: nrhs %" PetscInt_FMT ", sym %d, herm %d, size schur %" PetscInt_FMT ", size mat %" PetscInt_FMT "\n", solver, nrhs, symm, herm, size_schur, m));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* Test LU/Cholesky Factorization */
120c4762a1bSJed Brown   use_lu = PETSC_FALSE;
121c4762a1bSJed Brown   if (!symm) use_lu = PETSC_TRUE;
122c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
123c4762a1bSJed Brown   if (isolver == 1) use_lu = PETSC_TRUE;
124c4762a1bSJed Brown #endif
125c4762a1bSJed Brown   if (cuda && symm && !herm) use_lu = PETSC_TRUE;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   if (herm && !use_lu) { /* test also conversion routines inside the solver packages */
1289566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
1299566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATSEQSBAIJ, MAT_INPLACE_MATRIX, &A));
130c4762a1bSJed Brown   }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   if (use_lu) {
1339566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A, solver, MAT_FACTOR_LU, &F));
134c4762a1bSJed Brown   } else {
135c4762a1bSJed Brown     if (herm) {
1369566063dSJacob Faibussowitsch       PetscCall(MatSetOption(A, MAT_SPD, PETSC_TRUE));
137c4762a1bSJed Brown     } else {
1389566063dSJacob Faibussowitsch       PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
1399566063dSJacob Faibussowitsch       PetscCall(MatSetOption(A, MAT_SPD, PETSC_FALSE));
140c4762a1bSJed Brown     }
1419566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A, solver, MAT_FACTOR_CHOLESKY, &F));
142c4762a1bSJed Brown   }
1439566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, size_schur, m - size_schur, 1, &is_schur));
1449566063dSJacob Faibussowitsch   PetscCall(MatFactorSetSchurIS(F, is_schur));
145c4762a1bSJed Brown 
1469566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_schur));
147c4762a1bSJed Brown   if (use_lu) {
1489566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL));
149c4762a1bSJed Brown   } else {
1509566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL));
151c4762a1bSJed Brown   }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   for (nfact = 0; nfact < 3; nfact++) {
154c4762a1bSJed Brown     Mat AD;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown     if (!nfact) {
1579566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(x, rand));
158*48a46eb9SPierre Jolivet       if (symm && herm) PetscCall(VecAbs(x));
1599566063dSJacob Faibussowitsch       PetscCall(MatDiagonalSet(A, x, ADD_VALUES));
160c4762a1bSJed Brown     }
161c4762a1bSJed Brown     if (use_lu) {
1629566063dSJacob Faibussowitsch       PetscCall(MatLUFactorNumeric(F, A, NULL));
163c4762a1bSJed Brown     } else {
1649566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorNumeric(F, A, NULL));
165c4762a1bSJed Brown     }
166c4762a1bSJed Brown     if (cuda) {
1679566063dSJacob Faibussowitsch       PetscCall(MatFactorGetSchurComplement(F, &S, NULL));
1689566063dSJacob Faibussowitsch       PetscCall(MatSetType(S, MATSEQDENSECUDA));
1699566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(S, &xschur, &bschur));
1709566063dSJacob Faibussowitsch       PetscCall(MatFactorRestoreSchurComplement(F, &S, MAT_FACTOR_SCHUR_UNFACTORED));
171c4762a1bSJed Brown     }
1729566063dSJacob Faibussowitsch     PetscCall(MatFactorCreateSchurComplement(F, &S, NULL));
173*48a46eb9SPierre Jolivet     if (!cuda) PetscCall(MatCreateVecs(S, &xschur, &bschur));
1749566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xschur, &uschur));
175*48a46eb9SPierre Jolivet     if (nfact == 1 && (!cuda || (herm && symm))) PetscCall(MatFactorInvertSchurComplement(F));
176c4762a1bSJed Brown     for (nsolve = 0; nsolve < 2; nsolve++) {
1779566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(x, rand));
1789566063dSJacob Faibussowitsch       PetscCall(VecCopy(x, u));
179c4762a1bSJed Brown 
180c4762a1bSJed Brown       if (nsolve) {
1819566063dSJacob Faibussowitsch         PetscCall(MatMult(A, x, b));
1829566063dSJacob Faibussowitsch         PetscCall(MatSolve(F, b, x));
183c4762a1bSJed Brown       } else {
1849566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(A, x, b));
1859566063dSJacob Faibussowitsch         PetscCall(MatSolveTranspose(F, b, x));
186c4762a1bSJed Brown       }
187c4762a1bSJed Brown       /* Check the error */
1889566063dSJacob Faibussowitsch       PetscCall(VecAXPY(u, -1.0, x)); /* u <- (-1.0)x + u */
1899566063dSJacob Faibussowitsch       PetscCall(VecNorm(u, NORM_2, &norm));
190c4762a1bSJed Brown       if (norm > tol) {
191c4762a1bSJed Brown         PetscReal resi;
192c4762a1bSJed Brown         if (nsolve) {
1939566063dSJacob Faibussowitsch           PetscCall(MatMult(A, x, u)); /* u = A*x */
194c4762a1bSJed Brown         } else {
1959566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(A, x, u)); /* u = A*x */
196c4762a1bSJed Brown         }
1979566063dSJacob Faibussowitsch         PetscCall(VecAXPY(u, -1.0, b)); /* u <- (-1.0)b + u */
1989566063dSJacob Faibussowitsch         PetscCall(VecNorm(u, NORM_2, &resi));
199c4762a1bSJed Brown         if (nsolve) {
2009566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatSolve error: Norm of error %g, residual %g\n", nfact, nsolve, (double)norm, (double)resi));
201c4762a1bSJed Brown         } else {
2029566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatSolveTranspose error: Norm of error %g, residual %f\n", nfact, nsolve, (double)norm, (double)resi));
203c4762a1bSJed Brown         }
204c4762a1bSJed Brown       }
2059566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(xschur, rand));
2069566063dSJacob Faibussowitsch       PetscCall(VecCopy(xschur, uschur));
207c4762a1bSJed Brown       if (nsolve) {
2089566063dSJacob Faibussowitsch         PetscCall(MatMult(S, xschur, bschur));
2099566063dSJacob Faibussowitsch         PetscCall(MatFactorSolveSchurComplement(F, bschur, xschur));
210c4762a1bSJed Brown       } else {
2119566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(S, xschur, bschur));
2129566063dSJacob Faibussowitsch         PetscCall(MatFactorSolveSchurComplementTranspose(F, bschur, xschur));
213c4762a1bSJed Brown       }
214c4762a1bSJed Brown       /* Check the error */
2159566063dSJacob Faibussowitsch       PetscCall(VecAXPY(uschur, -1.0, xschur)); /* u <- (-1.0)x + u */
2169566063dSJacob Faibussowitsch       PetscCall(VecNorm(uschur, NORM_2, &norm));
217c4762a1bSJed Brown       if (norm > tol) {
218c4762a1bSJed Brown         PetscReal resi;
219c4762a1bSJed Brown         if (nsolve) {
2209566063dSJacob Faibussowitsch           PetscCall(MatMult(S, xschur, uschur)); /* u = A*x */
221c4762a1bSJed Brown         } else {
2229566063dSJacob Faibussowitsch           PetscCall(MatMultTranspose(S, xschur, uschur)); /* u = A*x */
223c4762a1bSJed Brown         }
2249566063dSJacob Faibussowitsch         PetscCall(VecAXPY(uschur, -1.0, bschur)); /* u <- (-1.0)b + u */
2259566063dSJacob Faibussowitsch         PetscCall(VecNorm(uschur, NORM_2, &resi));
226c4762a1bSJed Brown         if (nsolve) {
2279566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatFactorSolveSchurComplement error: Norm of error %g, residual %g\n", nfact, nsolve, (double)norm, (double)resi));
228c4762a1bSJed Brown         } else {
2299566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatFactorSolveSchurComplementTranspose error: Norm of error %g, residual %f\n", nfact, nsolve, (double)norm, (double)resi));
230c4762a1bSJed Brown         }
231c4762a1bSJed Brown       }
232c4762a1bSJed Brown     }
2339566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &AD));
234c4762a1bSJed Brown     if (!nfact) {
2359566063dSJacob Faibussowitsch       PetscCall(MatMatMult(AD, C, MAT_INITIAL_MATRIX, 2.0, &RHS));
236c4762a1bSJed Brown     } else {
2379566063dSJacob Faibussowitsch       PetscCall(MatMatMult(AD, C, MAT_REUSE_MATRIX, 2.0, &RHS));
238c4762a1bSJed Brown     }
2399566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AD));
240c4762a1bSJed Brown     for (nsolve = 0; nsolve < 2; nsolve++) {
2419566063dSJacob Faibussowitsch       PetscCall(MatMatSolve(F, RHS, X));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown       /* Check the error */
2449566063dSJacob Faibussowitsch       PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
2459566063dSJacob Faibussowitsch       PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
246*48a46eb9SPierre Jolivet       if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") MatMatSolve: Norm of error %g\n", nfact, nsolve, (double)norm));
247c4762a1bSJed Brown     }
248c4762a1bSJed Brown     if (isolver == 0) {
249c4762a1bSJed Brown       Mat spRHS, spRHST, RHST;
250c4762a1bSJed Brown 
2519566063dSJacob Faibussowitsch       PetscCall(MatTranspose(RHS, MAT_INITIAL_MATRIX, &RHST));
2529566063dSJacob Faibussowitsch       PetscCall(MatConvert(RHST, MATSEQAIJ, MAT_INITIAL_MATRIX, &spRHST));
2539566063dSJacob Faibussowitsch       PetscCall(MatCreateTranspose(spRHST, &spRHS));
254c4762a1bSJed Brown       for (nsolve = 0; nsolve < 2; nsolve++) {
2559566063dSJacob Faibussowitsch         PetscCall(MatMatSolve(F, spRHS, X));
256c4762a1bSJed Brown 
257c4762a1bSJed Brown         /* Check the error */
2589566063dSJacob Faibussowitsch         PetscCall(MatAXPY(X, -1.0, C, SAME_NONZERO_PATTERN));
2599566063dSJacob Faibussowitsch         PetscCall(MatNorm(X, NORM_FROBENIUS, &norm));
260*48a46eb9SPierre Jolivet         if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "(f %" PetscInt_FMT ", s %" PetscInt_FMT ") sparse MatMatSolve: Norm of error %g\n", nfact, nsolve, (double)norm));
261c4762a1bSJed Brown       }
2629566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&spRHST));
2639566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&spRHS));
2649566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&RHST));
265c4762a1bSJed Brown     }
2669566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&S));
2679566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&xschur));
2689566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bschur));
2699566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uschur));
270c4762a1bSJed Brown   }
271c4762a1bSJed Brown   /* Free data structures */
2729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
2749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
2759566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&X));
2769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&RHS));
2779566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
2789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
2809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2819566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
282b122ec5aSJacob Faibussowitsch   return 0;
283c4762a1bSJed Brown }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown /*TEST
286c4762a1bSJed Brown 
287c4762a1bSJed Brown    testset:
28838cf7977SStefano Zampini      requires: mkl_pardiso double !complex
289c4762a1bSJed Brown      args: -solver 1
290c4762a1bSJed Brown 
291c4762a1bSJed Brown      test:
292c4762a1bSJed Brown        suffix: mkl_pardiso
293c4762a1bSJed Brown      test:
294c4762a1bSJed Brown        requires: cuda
295c4762a1bSJed Brown        suffix: mkl_pardiso_cuda
296c4762a1bSJed Brown        args: -cuda_solve
297c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso.out
298c4762a1bSJed Brown      test:
299c4762a1bSJed Brown        suffix: mkl_pardiso_1
300c4762a1bSJed Brown        args: -symmetric_solve
301c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso_1.out
302c4762a1bSJed Brown      test:
303c4762a1bSJed Brown        requires: cuda
304c4762a1bSJed Brown        suffix: mkl_pardiso_cuda_1
305c4762a1bSJed Brown        args: -symmetric_solve -cuda_solve
306c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso_1.out
307c4762a1bSJed Brown      test:
308c4762a1bSJed Brown        suffix: mkl_pardiso_3
309c4762a1bSJed Brown        args: -symmetric_solve -hermitian_solve
310c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso_3.out
311c4762a1bSJed Brown      test:
312dfd57a17SPierre Jolivet        requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
313c4762a1bSJed Brown        suffix: mkl_pardiso_cuda_3
314c4762a1bSJed Brown        args: -symmetric_solve -hermitian_solve -cuda_solve
315c4762a1bSJed Brown        output_file: output/ex192_mkl_pardiso_3.out
316c4762a1bSJed Brown 
317c4762a1bSJed Brown    testset:
318c4762a1bSJed Brown      requires: mumps double !complex
319c4762a1bSJed Brown      args: -solver 0
320c4762a1bSJed Brown 
321c4762a1bSJed Brown      test:
322c4762a1bSJed Brown        suffix: mumps
323c4762a1bSJed Brown      test:
324c4762a1bSJed Brown        requires: cuda
325c4762a1bSJed Brown        suffix: mumps_cuda
326c4762a1bSJed Brown        args: -cuda_solve
327c4762a1bSJed Brown        output_file: output/ex192_mumps.out
328c4762a1bSJed Brown      test:
329c4762a1bSJed Brown        suffix: mumps_2
330c4762a1bSJed Brown        args: -symmetric_solve
331c4762a1bSJed Brown        output_file: output/ex192_mumps_2.out
332c4762a1bSJed Brown      test:
333c4762a1bSJed Brown        requires: cuda
334c4762a1bSJed Brown        suffix: mumps_cuda_2
335c4762a1bSJed Brown        args: -symmetric_solve -cuda_solve
336c4762a1bSJed Brown        output_file: output/ex192_mumps_2.out
337c4762a1bSJed Brown      test:
338c4762a1bSJed Brown        suffix: mumps_3
339c4762a1bSJed Brown        args: -symmetric_solve -hermitian_solve
340c4762a1bSJed Brown        output_file: output/ex192_mumps_3.out
341c4762a1bSJed Brown      test:
342dfd57a17SPierre Jolivet        requires: cuda defined(PETSC_HAVE_CUSOLVERDNDPOTRI)
343c4762a1bSJed Brown        suffix: mumps_cuda_3
344c4762a1bSJed Brown        args: -symmetric_solve -hermitian_solve -cuda_solve
345c4762a1bSJed Brown        output_file: output/ex192_mumps_3.out
346c4762a1bSJed Brown 
347c4762a1bSJed Brown TEST*/
348