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