1*18839329SStefano Zampini static char help[] = "Tests MatMatSolve() in Schur complement mode.\n\n"; 2*18839329SStefano Zampini 3*18839329SStefano Zampini #include <petscmat.h> 4*18839329SStefano Zampini 5*18839329SStefano Zampini int main(int argc, char **args) 6*18839329SStefano Zampini { 7*18839329SStefano Zampini Mat F, A, B, X, Y, S; 8*18839329SStefano Zampini IS is_schur; 9*18839329SStefano Zampini PetscMPIInt size; 10*18839329SStefano Zampini PetscInt ns = 0, m, n; 11*18839329SStefano Zampini PetscReal norm, tol = PETSC_SQRT_MACHINE_EPSILON; 12*18839329SStefano Zampini MatFactorType factor = MAT_FACTOR_LU; 13*18839329SStefano Zampini PetscViewer fd; 14*18839329SStefano Zampini char solver[256], converttype[256]; 15*18839329SStefano Zampini char file[PETSC_MAX_PATH_LEN]; 16*18839329SStefano Zampini PetscBool flg; 17*18839329SStefano Zampini 18*18839329SStefano Zampini PetscFunctionBeginUser; 19*18839329SStefano Zampini PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 20*18839329SStefano Zampini PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 21*18839329SStefano Zampini PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor test"); 22*18839329SStefano Zampini 23*18839329SStefano Zampini PetscCall(PetscOptionsGetString(NULL, NULL, "-A", file, sizeof(file), &flg)); 24*18839329SStefano Zampini PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Must provide a binary matrix with -A filename option"); 25*18839329SStefano Zampini PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 26*18839329SStefano Zampini PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 27*18839329SStefano Zampini PetscCall(MatLoad(A, fd)); 28*18839329SStefano Zampini PetscCall(PetscViewerDestroy(&fd)); 29*18839329SStefano Zampini PetscCall(MatGetSize(A, &m, &n)); 30*18839329SStefano Zampini PetscCheck(m == n, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n); 31*18839329SStefano Zampini PetscCall(MatViewFromOptions(A, NULL, "-A_view")); 32*18839329SStefano Zampini 33*18839329SStefano Zampini PetscCall(PetscOptionsGetString(NULL, NULL, "-B", file, sizeof(file), &flg)); 34*18839329SStefano Zampini PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Must provide a binary matrix with -B filename option"); 35*18839329SStefano Zampini PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 36*18839329SStefano Zampini PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 37*18839329SStefano Zampini PetscCall(MatLoad(B, fd)); 38*18839329SStefano Zampini PetscCall(PetscViewerDestroy(&fd)); 39*18839329SStefano Zampini PetscCall(MatViewFromOptions(B, NULL, "-B_view")); 40*18839329SStefano Zampini PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 41*18839329SStefano Zampini if (!flg) PetscCall(PetscObjectTypeCompare((PetscObject)B, MATTRANSPOSEVIRTUAL, &flg)); 42*18839329SStefano Zampini if (!flg) { 43*18839329SStefano Zampini Mat Bt; 44*18839329SStefano Zampini 45*18839329SStefano Zampini PetscCall(MatCreateTranspose(B, &Bt)); 46*18839329SStefano Zampini PetscCall(MatDestroy(&B)); 47*18839329SStefano Zampini B = Bt; 48*18839329SStefano Zampini } 49*18839329SStefano Zampini PetscCall(PetscOptionsGetString(NULL, NULL, "-B_convert_type", converttype, sizeof(converttype), &flg)); 50*18839329SStefano Zampini if (flg) PetscCall(MatConvert(B, converttype, MAT_INPLACE_MATRIX, &B)); 51*18839329SStefano Zampini 52*18839329SStefano Zampini PetscCall(PetscOptionsGetInt(NULL, NULL, "-ns", &ns, NULL)); 53*18839329SStefano Zampini 54*18839329SStefano Zampini PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_solver_type", solver, sizeof(solver), &flg)); 55*18839329SStefano Zampini if (!flg) PetscCall(PetscStrncpy(solver, MATSOLVERMUMPS, sizeof(solver))); 56*18839329SStefano Zampini PetscCall(PetscOptionsGetEnum(NULL, NULL, "-mat_factor_type", MatFactorTypes, (PetscEnum *)&factor, NULL)); 57*18839329SStefano Zampini PetscCall(MatGetFactor(A, solver, factor, &F)); 58*18839329SStefano Zampini 59*18839329SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF, ns, m - ns, 1, &is_schur)); 60*18839329SStefano Zampini PetscCall(MatFactorSetSchurIS(F, is_schur)); 61*18839329SStefano Zampini PetscCall(ISDestroy(&is_schur)); 62*18839329SStefano Zampini switch (factor) { 63*18839329SStefano Zampini case MAT_FACTOR_LU: 64*18839329SStefano Zampini PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL)); 65*18839329SStefano Zampini PetscCall(MatLUFactorNumeric(F, A, NULL)); 66*18839329SStefano Zampini break; 67*18839329SStefano Zampini case MAT_FACTOR_CHOLESKY: 68*18839329SStefano Zampini PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL)); 69*18839329SStefano Zampini PetscCall(MatCholeskyFactorNumeric(F, A, NULL)); 70*18839329SStefano Zampini break; 71*18839329SStefano Zampini default: 72*18839329SStefano Zampini PetscCheck(PETSC_FALSE, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not coded for factor type %s", MatFactorTypes[factor]); 73*18839329SStefano Zampini } 74*18839329SStefano Zampini 75*18839329SStefano Zampini PetscCall(MatFactorCreateSchurComplement(F, &S, NULL)); 76*18839329SStefano Zampini PetscCall(MatViewFromOptions(S, NULL, "-S_view")); 77*18839329SStefano Zampini PetscCall(MatDestroy(&S)); 78*18839329SStefano Zampini 79*18839329SStefano Zampini PetscCall(MatGetSize(B, NULL, &n)); 80*18839329SStefano Zampini PetscCall(MatCreate(PETSC_COMM_WORLD, &X)); 81*18839329SStefano Zampini PetscCall(MatSetSizes(X, m, PETSC_DECIDE, PETSC_DECIDE, n)); 82*18839329SStefano Zampini PetscCall(MatSetType(X, MATDENSE)); 83*18839329SStefano Zampini PetscCall(MatSetFromOptions(X)); 84*18839329SStefano Zampini PetscCall(MatSetUp(X)); 85*18839329SStefano Zampini 86*18839329SStefano Zampini PetscCall(MatMatSolve(F, B, X)); 87*18839329SStefano Zampini PetscCall(MatViewFromOptions(X, NULL, "-X_view")); 88*18839329SStefano Zampini PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Y)); 89*18839329SStefano Zampini PetscCall(MatViewFromOptions(Y, NULL, "-Y_view")); 90*18839329SStefano Zampini PetscCall(MatAXPY(Y, -1.0, B, SAME_NONZERO_PATTERN)); 91*18839329SStefano Zampini PetscCall(MatViewFromOptions(Y, NULL, "-err_view")); 92*18839329SStefano Zampini PetscCall(MatNorm(Y, NORM_FROBENIUS, &norm)); 93*18839329SStefano Zampini if (norm > tol) { 94*18839329SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMatSolve: Norm of error %g\n", (double)norm)); 95*18839329SStefano Zampini PetscCall(MatConvert(Y, MATAIJ, MAT_INPLACE_MATRIX, &Y)); 96*18839329SStefano Zampini PetscCall(MatFilter(Y, PETSC_SMALL, PETSC_TRUE, PETSC_FALSE)); 97*18839329SStefano Zampini PetscCall(MatViewFromOptions(Y, NULL, "-aij_err_view")); 98*18839329SStefano Zampini } 99*18839329SStefano Zampini PetscCall(MatDestroy(&A)); 100*18839329SStefano Zampini PetscCall(MatDestroy(&X)); 101*18839329SStefano Zampini PetscCall(MatDestroy(&F)); 102*18839329SStefano Zampini PetscCall(MatDestroy(&B)); 103*18839329SStefano Zampini PetscCall(MatDestroy(&Y)); 104*18839329SStefano Zampini PetscCall(PetscFinalize()); 105*18839329SStefano Zampini return 0; 106*18839329SStefano Zampini } 107*18839329SStefano Zampini 108*18839329SStefano Zampini /*TEST 109*18839329SStefano Zampini 110*18839329SStefano Zampini test: 111*18839329SStefano Zampini output_file: output/ex62_1.out 112*18839329SStefano Zampini suffix: mumps_1 113*18839329SStefano Zampini requires: mumps double !complex !defined(PETSC_USE_64BIT_INDICES) 114*18839329SStefano Zampini args: -A ${DATAFILESPATH}/matrices/factorSchur/A.dat -B ${DATAFILESPATH}/matrices/factorSchur/B1.dat -ns {{0 1}} 115*18839329SStefano Zampini 116*18839329SStefano Zampini test: 117*18839329SStefano Zampini output_file: output/ex62_1.out 118*18839329SStefano Zampini suffix: mumps_2 119*18839329SStefano Zampini requires: mumps double !complex !defined(PETSC_USE_64BIT_INDICES) 120*18839329SStefano Zampini args: -A ${DATAFILESPATH}/matrices/factorSchur/A.dat -B ${DATAFILESPATH}/matrices/factorSchur/B2.dat -ns {{0 1}} 121*18839329SStefano Zampini 122*18839329SStefano Zampini TEST*/ 123