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