xref: /petsc/src/mat/tests/ex267.c (revision ae1ee55146a7ad071171b861759b85940c7e5c67)
1f5bab676SJose E. Roman static char help[] = "Test different MatSolve routines with MATTRANSPOSEVIRTUAL.\n\n";
2f5bab676SJose E. Roman 
3f5bab676SJose E. Roman #include <petscmat.h>
4f5bab676SJose E. Roman 
TestMatrix(const char * test,Mat A,PetscInt nrhs,PetscBool inplace,PetscBool chol)5f5bab676SJose E. Roman PetscErrorCode TestMatrix(const char *test, Mat A, PetscInt nrhs, PetscBool inplace, PetscBool chol)
6f5bab676SJose E. Roman {
7f5bab676SJose E. Roman   Mat       F, RHS, X, C1;
8f5bab676SJose E. Roman   Vec       b, x, y, f;
9f5bab676SJose E. Roman   IS        perm, iperm;
10f5bab676SJose E. Roman   PetscInt  n, i;
11f5bab676SJose E. Roman   PetscReal norm, tol = 1000 * PETSC_MACHINE_EPSILON;
12f5bab676SJose E. Roman   PetscBool ht;
13f5bab676SJose E. Roman #if defined(PETSC_USE_COMPLEX)
14f5bab676SJose E. Roman   PetscScalar v1 = PetscCMPLX(1.0, -0.1), v2 = PetscCMPLX(-1.0, 0.1);
15f5bab676SJose E. Roman #else
16f5bab676SJose E. Roman   PetscScalar v1 = 1.0, v2 = -1.0;
17f5bab676SJose E. Roman #endif
18f5bab676SJose E. Roman 
19f5bab676SJose E. Roman   PetscFunctionBegin;
20f5bab676SJose E. Roman   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATHERMITIANTRANSPOSEVIRTUAL, &ht));
21f5bab676SJose E. Roman   PetscCall(MatCreateVecs(A, &f, &b));
22f5bab676SJose E. Roman   PetscCall(MatCreateVecs(A, &x, &y));
23f5bab676SJose E. Roman   PetscCall(VecSet(b, v1));
24f5bab676SJose E. Roman   PetscCall(VecSet(y, v2));
25f5bab676SJose E. Roman 
26f5bab676SJose E. Roman   PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm));
27f5bab676SJose E. Roman   if (!inplace) {
28f5bab676SJose E. Roman     if (!chol) {
29f5bab676SJose E. Roman       PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &F));
30f5bab676SJose E. Roman       PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, NULL));
31f5bab676SJose E. Roman       PetscCall(MatLUFactorNumeric(F, A, NULL));
32f5bab676SJose E. Roman     } else { /* Cholesky */
33f5bab676SJose E. Roman       PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &F));
34f5bab676SJose E. Roman       PetscCall(MatCholeskyFactorSymbolic(F, A, perm, NULL));
35f5bab676SJose E. Roman       PetscCall(MatCholeskyFactorNumeric(F, A, NULL));
36f5bab676SJose E. Roman     }
37f5bab676SJose E. Roman   } else { /* Test inplace factorization */
38f5bab676SJose E. Roman     PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &F));
39f5bab676SJose E. Roman     if (!chol) PetscCall(MatLUFactor(F, perm, iperm, NULL));
40f5bab676SJose E. Roman     else PetscCall(MatCholeskyFactor(F, perm, NULL));
41f5bab676SJose E. Roman   }
42f5bab676SJose E. Roman 
43f5bab676SJose E. Roman   /* MatSolve */
44f5bab676SJose E. Roman   PetscCall(MatSolve(F, b, x));
45f5bab676SJose E. Roman   PetscCall(MatMult(A, x, f));
46f5bab676SJose E. Roman   PetscCall(VecAXPY(f, -1.0, b));
47f5bab676SJose E. Roman   PetscCall(VecNorm(f, NORM_2, &norm));
48f5bab676SJose E. Roman   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%12s MatSolve              : Error of norm %g\n", test, (double)norm));
49f5bab676SJose E. Roman 
50f5bab676SJose E. Roman   /* MatSolveTranspose */
51f5bab676SJose E. Roman   if (!ht) {
52f5bab676SJose E. Roman     PetscCall(MatSolveTranspose(F, b, x));
53f5bab676SJose E. Roman     PetscCall(MatMultTranspose(A, x, f));
54f5bab676SJose E. Roman     PetscCall(VecAXPY(f, -1.0, b));
55f5bab676SJose E. Roman     PetscCall(VecNorm(f, NORM_2, &norm));
56f5bab676SJose E. Roman     if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%12s MatSolveTranspose     : Error of norm %g\n", test, (double)norm));
57f5bab676SJose E. Roman   }
58f5bab676SJose E. Roman 
59f5bab676SJose E. Roman   /* MatSolveAdd */
60f5bab676SJose E. Roman   PetscCall(MatSolveAdd(F, b, y, x));
61f5bab676SJose E. Roman   PetscCall(MatMult(A, y, f));
62f5bab676SJose E. Roman   PetscCall(VecScale(f, -1.0));
63f5bab676SJose E. Roman   PetscCall(MatMultAdd(A, x, f, f));
64f5bab676SJose E. Roman   PetscCall(VecAXPY(f, -1.0, b));
65f5bab676SJose E. Roman   PetscCall(VecNorm(f, NORM_2, &norm));
66f5bab676SJose E. Roman   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%12s MatSolveAdd           : Error of norm %g\n", test, (double)norm));
67f5bab676SJose E. Roman 
68f5bab676SJose E. Roman   /* MatSolveTransposeAdd */
69f5bab676SJose E. Roman   if (!ht) {
70f5bab676SJose E. Roman     PetscCall(MatSolveTransposeAdd(F, b, y, x));
71f5bab676SJose E. Roman     PetscCall(MatMultTranspose(A, y, f));
72f5bab676SJose E. Roman     PetscCall(VecScale(f, -1.0));
73f5bab676SJose E. Roman     PetscCall(MatMultTransposeAdd(A, x, f, f));
74f5bab676SJose E. Roman     PetscCall(VecAXPY(f, -1.0, b));
75f5bab676SJose E. Roman     PetscCall(VecNorm(f, NORM_2, &norm));
76f5bab676SJose E. Roman     if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%12s MatSolveTransposeAdd  : Error of norm %g\n", test, (double)norm));
77f5bab676SJose E. Roman   }
78f5bab676SJose E. Roman 
79f5bab676SJose E. Roman   /* MatMatSolve */
80f5bab676SJose E. Roman   PetscCall(MatGetSize(A, &n, NULL));
81f5bab676SJose E. Roman   PetscCall(MatCreate(PETSC_COMM_WORLD, &RHS));
82f5bab676SJose E. Roman   PetscCall(MatSetSizes(RHS, PETSC_DECIDE, PETSC_DECIDE, n, nrhs));
83f5bab676SJose E. Roman   PetscCall(MatSetType(RHS, MATSEQDENSE));
84f5bab676SJose E. Roman   PetscCall(MatSetUp(RHS));
85f5bab676SJose E. Roman   for (i = 0; i < nrhs; i++) PetscCall(MatSetValue(RHS, i, i, 1.0, INSERT_VALUES));
86f5bab676SJose E. Roman   PetscCall(MatAssemblyBegin(RHS, MAT_FINAL_ASSEMBLY));
87f5bab676SJose E. Roman   PetscCall(MatAssemblyEnd(RHS, MAT_FINAL_ASSEMBLY));
88f5bab676SJose E. Roman   PetscCall(MatDuplicate(RHS, MAT_DO_NOT_COPY_VALUES, &X));
89f5bab676SJose E. Roman 
90f5bab676SJose E. Roman   if (!ht) {
91f5bab676SJose E. Roman     PetscCall(MatMatSolve(F, RHS, X));
92fb842aefSJose E. Roman     PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &C1));
93f5bab676SJose E. Roman     PetscCall(MatAXPY(C1, -1.0, RHS, SAME_NONZERO_PATTERN));
94f5bab676SJose E. Roman     PetscCall(MatNorm(C1, NORM_FROBENIUS, &norm));
95f5bab676SJose E. Roman     if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%12s MatMatSolve           : Error of norm %g\n", test, (double)norm));
96f5bab676SJose E. Roman     PetscCall(MatDestroy(&C1));
97f5bab676SJose E. Roman 
98f5bab676SJose E. Roman     PetscCall(MatMatSolveTranspose(F, RHS, X));
99fb842aefSJose E. Roman     PetscCall(MatTransposeMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &C1));
100f5bab676SJose E. Roman     PetscCall(MatAXPY(C1, -1.0, RHS, SAME_NONZERO_PATTERN));
101f5bab676SJose E. Roman     PetscCall(MatNorm(C1, NORM_FROBENIUS, &norm));
102f5bab676SJose E. Roman     if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%12s MatMatSolveTranspose  : Error of norm %g\n", test, (double)norm));
103f5bab676SJose E. Roman     PetscCall(MatDestroy(&C1));
104f5bab676SJose E. Roman   }
105f5bab676SJose E. Roman   PetscCall(VecDestroy(&b));
106f5bab676SJose E. Roman   PetscCall(VecDestroy(&x));
107f5bab676SJose E. Roman   PetscCall(VecDestroy(&f));
108f5bab676SJose E. Roman   PetscCall(VecDestroy(&y));
109f5bab676SJose E. Roman   PetscCall(ISDestroy(&perm));
110f5bab676SJose E. Roman   PetscCall(ISDestroy(&iperm));
111f5bab676SJose E. Roman   PetscCall(MatDestroy(&F));
112f5bab676SJose E. Roman   PetscCall(MatDestroy(&RHS));
113f5bab676SJose E. Roman   PetscCall(MatDestroy(&X));
114f5bab676SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
115f5bab676SJose E. Roman }
116f5bab676SJose E. Roman 
main(int argc,char ** args)117f5bab676SJose E. Roman int main(int argc, char **args)
118f5bab676SJose E. Roman {
119f5bab676SJose E. Roman   PetscMPIInt size;
120f5bab676SJose E. Roman   Mat         A, At, Aht;
121f5bab676SJose E. Roman   PetscInt    i, n = 8, nrhs = 2;
122f5bab676SJose E. Roman   PetscBool   aij, inplace = PETSC_FALSE;
123f5bab676SJose E. Roman #if defined(PETSC_USE_COMPLEX)
124f5bab676SJose E. Roman   PetscScalar a = PetscCMPLX(-1.0, 0.5);
125f5bab676SJose E. Roman #else
126f5bab676SJose E. Roman   PetscScalar a = -1.0;
127f5bab676SJose E. Roman #endif
128f5bab676SJose E. Roman 
129f5bab676SJose E. Roman   PetscFunctionBeginUser;
130c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
131f5bab676SJose E. Roman   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
132f5bab676SJose E. Roman   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only");
133f5bab676SJose E. Roman   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
134f5bab676SJose E. Roman   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nrhs", &nrhs, NULL));
135f5bab676SJose E. Roman   PetscCall(PetscOptionsGetBool(NULL, NULL, "-inplace", &inplace, NULL));
136f5bab676SJose E. Roman   PetscCheck(nrhs <= n, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Must have nrhs <= n");
137f5bab676SJose E. Roman 
138f5bab676SJose E. Roman   /* Hermitian matrix */
139f5bab676SJose E. Roman   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
140f5bab676SJose E. Roman   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, n, n));
141f5bab676SJose E. Roman   PetscCall(MatSetFromOptions(A));
142f5bab676SJose E. Roman   for (i = 0; i < n; i++) {
143f5bab676SJose E. Roman     if (i > 0) PetscCall(MatSetValue(A, i, i - 1, a, INSERT_VALUES));
144f5bab676SJose E. Roman     if (i < n - 1) PetscCall(MatSetValue(A, i, i + 1, PetscConj(a), INSERT_VALUES));
145f5bab676SJose E. Roman     PetscCall(MatSetValue(A, i, i, 2.0, INSERT_VALUES));
146f5bab676SJose E. Roman   }
147f5bab676SJose E. Roman   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
148f5bab676SJose E. Roman   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
149f5bab676SJose E. Roman 
150f5bab676SJose E. Roman   PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &aij, MATSEQAIJ, MATSEQBAIJ, ""));
151*fc2fb351SPierre Jolivet   PetscCall(MatSetOption(A, PetscDefined(USE_COMPLEX) ? MAT_HERMITIAN : MAT_SYMMETRIC, PETSC_TRUE));
152f5bab676SJose E. Roman 
153f5bab676SJose E. Roman   PetscCall(MatCreateTranspose(A, &At));
154f5bab676SJose E. Roman   PetscCall(MatCreateHermitianTranspose(A, &Aht));
155f5bab676SJose E. Roman 
156f5bab676SJose E. Roman   PetscCall(TestMatrix("LU T", At, nrhs, inplace, PETSC_FALSE));
157f5bab676SJose E. Roman   PetscCall(TestMatrix("LU HT", Aht, nrhs, inplace, PETSC_FALSE));
158f5bab676SJose E. Roman   if (!aij) {
159f5bab676SJose E. Roman     PetscCall(TestMatrix("Chol T", At, nrhs, inplace, PETSC_TRUE));
160f5bab676SJose E. Roman     PetscCall(TestMatrix("Chol HT", Aht, nrhs, inplace, PETSC_TRUE));
161f5bab676SJose E. Roman   }
162f5bab676SJose E. Roman 
163f5bab676SJose E. Roman   /* Make the matrix non-Hermitian */
164f5bab676SJose E. Roman   PetscCall(MatSetValue(A, 0, 1, -5.0, INSERT_VALUES));
165f5bab676SJose E. Roman   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
166f5bab676SJose E. Roman   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
167f5bab676SJose E. Roman #if defined(PETSC_USE_COMPLEX)
168f5bab676SJose E. Roman   PetscCall(MatSetOption(A, MAT_HERMITIAN, PETSC_FALSE));
169f5bab676SJose E. Roman #else
170f5bab676SJose E. Roman   PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_FALSE));
171f5bab676SJose E. Roman #endif
172f5bab676SJose E. Roman 
173f5bab676SJose E. Roman   PetscCall(TestMatrix("LU T nonsym", At, nrhs, inplace, PETSC_FALSE));
174f5bab676SJose E. Roman   PetscCall(TestMatrix("LU HT nonsym", Aht, nrhs, inplace, PETSC_FALSE));
175f5bab676SJose E. Roman 
176f5bab676SJose E. Roman   PetscCall(MatDestroy(&A));
177f5bab676SJose E. Roman   PetscCall(MatDestroy(&At));
178f5bab676SJose E. Roman   PetscCall(MatDestroy(&Aht));
179f5bab676SJose E. Roman   PetscCall(PetscFinalize());
180f5bab676SJose E. Roman   return 0;
181f5bab676SJose E. Roman }
182f5bab676SJose E. Roman 
183f5bab676SJose E. Roman /*TEST
184f5bab676SJose E. Roman 
185f5bab676SJose E. Roman    test:
186f5bab676SJose E. Roman       suffix: 1
187f5bab676SJose E. Roman       args: -inplace {{0 1}} -mat_type {{aij dense}}
188f5bab676SJose E. Roman       output_file: output/empty.out
189f5bab676SJose E. Roman 
190f5bab676SJose E. Roman TEST*/
191