1bd5dbebeSNils Friess static char help[] = "Tests MatForwardSolve and MatBackwardSolve for LU and Cholesky decompositions using C/Pardiso.\n\n";
2bd5dbebeSNils Friess
3bd5dbebeSNils Friess #include <petscdmda.h>
4bd5dbebeSNils Friess #include <petscmat.h>
5bd5dbebeSNils Friess
main(int argc,char ** args)6bd5dbebeSNils Friess int main(int argc, char **args)
7bd5dbebeSNils Friess {
8bd5dbebeSNils Friess DM da;
9bd5dbebeSNils Friess DMDALocalInfo info;
10bd5dbebeSNils Friess Mat A, F;
11bd5dbebeSNils Friess Vec x, y, b, ytmp;
12bd5dbebeSNils Friess IS rowis, colis;
13bd5dbebeSNils Friess PetscInt i, j, k, n = 5;
14bd5dbebeSNils Friess PetscBool CHOL = PETSC_FALSE;
15bd5dbebeSNils Friess MatStencil row, cols[5];
16bd5dbebeSNils Friess PetscScalar vals[5];
17bd5dbebeSNils Friess PetscReal norm2, tol = 100. * PETSC_MACHINE_EPSILON;
18bd5dbebeSNils Friess PetscRandom rdm;
19bd5dbebeSNils Friess MatFactorInfo finfo;
20bd5dbebeSNils Friess
21bd5dbebeSNils Friess PetscFunctionBeginUser;
22*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
23bd5dbebeSNils Friess PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
24bd5dbebeSNils Friess
25bd5dbebeSNils Friess PetscCall(DMDACreate2d(MPI_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, n, n, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da));
26bd5dbebeSNils Friess PetscCall(DMSetUp(da));
27bd5dbebeSNils Friess
28bd5dbebeSNils Friess PetscCall(DMCreateMatrix(da, &A));
29bd5dbebeSNils Friess PetscCall(PetscOptionsHasName(NULL, NULL, "-chol", &CHOL));
30bd5dbebeSNils Friess if (CHOL) PetscCall(MatSetType(A, MATSBAIJ));
31bd5dbebeSNils Friess else PetscCall(MatSetType(A, MATBAIJ));
32bd5dbebeSNils Friess PetscCall(MatSetUp(A));
33bd5dbebeSNils Friess
34bd5dbebeSNils Friess PetscCall(DMDAGetLocalInfo(da, &info));
35bd5dbebeSNils Friess for (j = info.ys; j < info.ys + info.ym; j++) {
36bd5dbebeSNils Friess for (i = info.xs; i < info.xs + info.xm; i++) {
37bd5dbebeSNils Friess row.j = j;
38bd5dbebeSNils Friess row.i = i;
39bd5dbebeSNils Friess
40bd5dbebeSNils Friess k = 0;
41bd5dbebeSNils Friess if (j != 0) {
42bd5dbebeSNils Friess cols[k].j = j - 1;
43bd5dbebeSNils Friess cols[k].i = i;
44bd5dbebeSNils Friess vals[k] = -1;
45bd5dbebeSNils Friess ++k;
46bd5dbebeSNils Friess }
47bd5dbebeSNils Friess if (i != 0) {
48bd5dbebeSNils Friess cols[k].j = j;
49bd5dbebeSNils Friess cols[k].i = i - 1;
50bd5dbebeSNils Friess vals[k] = -1;
51bd5dbebeSNils Friess ++k;
52bd5dbebeSNils Friess }
53bd5dbebeSNils Friess cols[k].j = j;
54bd5dbebeSNils Friess cols[k].i = i;
55bd5dbebeSNils Friess vals[k] = 4;
56bd5dbebeSNils Friess ++k;
57bd5dbebeSNils Friess if (j != info.my - 1) {
58bd5dbebeSNils Friess cols[k].j = j + 1;
59bd5dbebeSNils Friess cols[k].i = i;
60bd5dbebeSNils Friess vals[k] = -1;
61bd5dbebeSNils Friess ++k;
62bd5dbebeSNils Friess }
63bd5dbebeSNils Friess if (i != info.mx - 1) {
64bd5dbebeSNils Friess cols[k].j = j;
65bd5dbebeSNils Friess cols[k].i = i + 1;
66bd5dbebeSNils Friess vals[k] = -1;
67bd5dbebeSNils Friess ++k;
68bd5dbebeSNils Friess }
69bd5dbebeSNils Friess PetscCall(MatSetValuesStencil(A, 1, &row, k, cols, vals, INSERT_VALUES));
70bd5dbebeSNils Friess }
71bd5dbebeSNils Friess }
72bd5dbebeSNils Friess PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
73bd5dbebeSNils Friess PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
74bd5dbebeSNils Friess
75bd5dbebeSNils Friess if (CHOL) PetscCall(MatSetOption(A, MAT_SPD, PETSC_TRUE));
76bd5dbebeSNils Friess
77bd5dbebeSNils Friess /* Create vectors for error checking */
78bd5dbebeSNils Friess PetscCall(MatCreateVecs(A, &x, &b));
79bd5dbebeSNils Friess PetscCall(VecDuplicate(x, &y));
80bd5dbebeSNils Friess PetscCall(VecDuplicate(x, &ytmp));
81bd5dbebeSNils Friess PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
82bd5dbebeSNils Friess PetscCall(PetscRandomSetFromOptions(rdm));
83bd5dbebeSNils Friess PetscCall(VecSetRandom(x, rdm));
84bd5dbebeSNils Friess PetscCall(MatMult(A, x, b));
85bd5dbebeSNils Friess
86bd5dbebeSNils Friess PetscCall(MatGetOrdering(A, MATORDERINGNATURAL, &rowis, &colis));
87bd5dbebeSNils Friess if (CHOL) {
88bd5dbebeSNils Friess PetscCall(PetscPrintf(MPI_COMM_WORLD, "Test Cholesky...\n"));
89bd5dbebeSNils Friess PetscCall(MatGetFactor(A, MATSOLVERMKL_CPARDISO, MAT_FACTOR_CHOLESKY, &F));
90bd5dbebeSNils Friess PetscCall(MatCholeskyFactorSymbolic(F, A, rowis, &finfo));
91bd5dbebeSNils Friess PetscCall(MatCholeskyFactorNumeric(F, A, &finfo));
92bd5dbebeSNils Friess } else {
93bd5dbebeSNils Friess PetscCall(PetscPrintf(MPI_COMM_WORLD, "Test LU...\n"));
94bd5dbebeSNils Friess PetscCall(MatGetFactor(A, MATSOLVERMKL_CPARDISO, MAT_FACTOR_LU, &F));
95bd5dbebeSNils Friess PetscCall(MatLUFactorSymbolic(F, A, rowis, colis, &finfo));
96bd5dbebeSNils Friess PetscCall(MatLUFactorNumeric(F, A, &finfo));
97bd5dbebeSNils Friess }
98bd5dbebeSNils Friess
99bd5dbebeSNils Friess PetscCall(PetscPrintf(MPI_COMM_WORLD, "Test MatForwardSolve...\n"));
100bd5dbebeSNils Friess PetscCall(MatForwardSolve(F, b, ytmp));
101bd5dbebeSNils Friess PetscCall(PetscPrintf(MPI_COMM_WORLD, "Test MatBackwardSolve...\n"));
102bd5dbebeSNils Friess PetscCall(MatBackwardSolve(F, ytmp, y));
103bd5dbebeSNils Friess PetscCall(VecAXPY(y, -1.0, x));
104bd5dbebeSNils Friess PetscCall(VecNorm(y, NORM_2, &norm2));
105bd5dbebeSNils Friess if (norm2 > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatForwardSolve and BackwardSolve: Norm of error=%g\n", (double)norm2));
106bd5dbebeSNils Friess
107bd5dbebeSNils Friess /* Free data structures */
108bd5dbebeSNils Friess PetscCall(ISDestroy(&rowis));
109bd5dbebeSNils Friess PetscCall(ISDestroy(&colis));
110bd5dbebeSNils Friess PetscCall(MatDestroy(&F));
111bd5dbebeSNils Friess PetscCall(MatDestroy(&A));
112bd5dbebeSNils Friess PetscCall(DMDestroy(&da));
113bd5dbebeSNils Friess PetscCall(PetscRandomDestroy(&rdm));
114bd5dbebeSNils Friess PetscCall(VecDestroy(&x));
115bd5dbebeSNils Friess PetscCall(VecDestroy(&y));
116bd5dbebeSNils Friess PetscCall(VecDestroy(&ytmp));
117bd5dbebeSNils Friess PetscCall(VecDestroy(&b));
118bd5dbebeSNils Friess PetscCall(PetscFinalize());
119bd5dbebeSNils Friess return 0;
120bd5dbebeSNils Friess }
121bd5dbebeSNils Friess
122bd5dbebeSNils Friess /*TEST
123bd5dbebeSNils Friess
124bd5dbebeSNils Friess test:
125bd5dbebeSNils Friess requires: mkl_cpardiso
126bd5dbebeSNils Friess nsize: 4
127bd5dbebeSNils Friess
128bd5dbebeSNils Friess test:
129bd5dbebeSNils Friess suffix: 2
130bd5dbebeSNils Friess requires: !complex mkl_cpardiso
131bd5dbebeSNils Friess nsize: 4
132bd5dbebeSNils Friess args: -chol
133bd5dbebeSNils Friess
134bd5dbebeSNils Friess TEST*/
135