118839329SStefano Zampini static char help[] = "Tests MatMatSolve() in Schur complement mode.\n\n"; 218839329SStefano Zampini 318839329SStefano Zampini #include <petscmat.h> 418839329SStefano Zampini 518839329SStefano Zampini int main(int argc, char **args) 618839329SStefano Zampini { 718839329SStefano Zampini Mat F, A, B, X, Y, S; 818839329SStefano Zampini IS is_schur; 918839329SStefano Zampini PetscMPIInt size; 1018839329SStefano Zampini PetscInt ns = 0, m, n; 1118839329SStefano Zampini PetscReal norm, tol = PETSC_SQRT_MACHINE_EPSILON; 1218839329SStefano Zampini MatFactorType factor = MAT_FACTOR_LU; 1318839329SStefano Zampini PetscViewer fd; 1418839329SStefano Zampini char solver[256], converttype[256]; 1518839329SStefano Zampini char file[PETSC_MAX_PATH_LEN]; 1618839329SStefano Zampini PetscBool flg; 1718839329SStefano Zampini 1818839329SStefano Zampini PetscFunctionBeginUser; 1918839329SStefano Zampini PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 2018839329SStefano Zampini PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2118839329SStefano Zampini PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor test"); 2218839329SStefano Zampini 2318839329SStefano Zampini PetscCall(PetscOptionsGetString(NULL, NULL, "-A", file, sizeof(file), &flg)); 2418839329SStefano Zampini PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Must provide a binary matrix with -A filename option"); 2518839329SStefano Zampini PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 2618839329SStefano Zampini PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 2718839329SStefano Zampini PetscCall(MatLoad(A, fd)); 2818839329SStefano Zampini PetscCall(PetscViewerDestroy(&fd)); 2918839329SStefano Zampini PetscCall(MatGetSize(A, &m, &n)); 3018839329SStefano 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); 3118839329SStefano Zampini PetscCall(MatViewFromOptions(A, NULL, "-A_view")); 3218839329SStefano Zampini 3318839329SStefano Zampini PetscCall(PetscOptionsGetString(NULL, NULL, "-B", file, sizeof(file), &flg)); 3418839329SStefano Zampini PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Must provide a binary matrix with -B filename option"); 3518839329SStefano Zampini PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 3618839329SStefano Zampini PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 3718839329SStefano Zampini PetscCall(MatLoad(B, fd)); 3818839329SStefano Zampini PetscCall(PetscViewerDestroy(&fd)); 3918839329SStefano Zampini PetscCall(MatViewFromOptions(B, NULL, "-B_view")); 4018839329SStefano Zampini PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)B, &flg, MATSEQDENSE, MATMPIDENSE, NULL)); 4118839329SStefano Zampini if (!flg) PetscCall(PetscObjectTypeCompare((PetscObject)B, MATTRANSPOSEVIRTUAL, &flg)); 4218839329SStefano Zampini if (!flg) { 4318839329SStefano Zampini Mat Bt; 4418839329SStefano Zampini 4518839329SStefano Zampini PetscCall(MatCreateTranspose(B, &Bt)); 4618839329SStefano Zampini PetscCall(MatDestroy(&B)); 4718839329SStefano Zampini B = Bt; 4818839329SStefano Zampini } 4918839329SStefano Zampini PetscCall(PetscOptionsGetString(NULL, NULL, "-B_convert_type", converttype, sizeof(converttype), &flg)); 5018839329SStefano Zampini if (flg) PetscCall(MatConvert(B, converttype, MAT_INPLACE_MATRIX, &B)); 5118839329SStefano Zampini 5218839329SStefano Zampini PetscCall(PetscOptionsGetInt(NULL, NULL, "-ns", &ns, NULL)); 5318839329SStefano Zampini 5418839329SStefano Zampini PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_solver_type", solver, sizeof(solver), &flg)); 5518839329SStefano Zampini if (!flg) PetscCall(PetscStrncpy(solver, MATSOLVERMUMPS, sizeof(solver))); 5618839329SStefano Zampini PetscCall(PetscOptionsGetEnum(NULL, NULL, "-mat_factor_type", MatFactorTypes, (PetscEnum *)&factor, NULL)); 5718839329SStefano Zampini PetscCall(MatGetFactor(A, solver, factor, &F)); 5818839329SStefano Zampini 5918839329SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF, ns, m - ns, 1, &is_schur)); 6018839329SStefano Zampini PetscCall(MatFactorSetSchurIS(F, is_schur)); 6118839329SStefano Zampini PetscCall(ISDestroy(&is_schur)); 6218839329SStefano Zampini switch (factor) { 6318839329SStefano Zampini case MAT_FACTOR_LU: 6418839329SStefano Zampini PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL)); 6518839329SStefano Zampini PetscCall(MatLUFactorNumeric(F, A, NULL)); 6618839329SStefano Zampini break; 6718839329SStefano Zampini case MAT_FACTOR_CHOLESKY: 6818839329SStefano Zampini PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL)); 6918839329SStefano Zampini PetscCall(MatCholeskyFactorNumeric(F, A, NULL)); 7018839329SStefano Zampini break; 7118839329SStefano Zampini default: 7218839329SStefano Zampini PetscCheck(PETSC_FALSE, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not coded for factor type %s", MatFactorTypes[factor]); 7318839329SStefano Zampini } 7418839329SStefano Zampini 7518839329SStefano Zampini PetscCall(MatFactorCreateSchurComplement(F, &S, NULL)); 7618839329SStefano Zampini PetscCall(MatViewFromOptions(S, NULL, "-S_view")); 7718839329SStefano Zampini PetscCall(MatDestroy(&S)); 7818839329SStefano Zampini 7918839329SStefano Zampini PetscCall(MatGetSize(B, NULL, &n)); 8018839329SStefano Zampini PetscCall(MatCreate(PETSC_COMM_WORLD, &X)); 8118839329SStefano Zampini PetscCall(MatSetSizes(X, m, PETSC_DECIDE, PETSC_DECIDE, n)); 8218839329SStefano Zampini PetscCall(MatSetType(X, MATDENSE)); 8318839329SStefano Zampini PetscCall(MatSetFromOptions(X)); 8418839329SStefano Zampini PetscCall(MatSetUp(X)); 8518839329SStefano Zampini 8618839329SStefano Zampini PetscCall(MatMatSolve(F, B, X)); 8718839329SStefano Zampini PetscCall(MatViewFromOptions(X, NULL, "-X_view")); 88fb842aefSJose E. Roman PetscCall(MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Y)); 8918839329SStefano Zampini PetscCall(MatViewFromOptions(Y, NULL, "-Y_view")); 9018839329SStefano Zampini PetscCall(MatAXPY(Y, -1.0, B, SAME_NONZERO_PATTERN)); 9118839329SStefano Zampini PetscCall(MatViewFromOptions(Y, NULL, "-err_view")); 9218839329SStefano Zampini PetscCall(MatNorm(Y, NORM_FROBENIUS, &norm)); 9318839329SStefano Zampini if (norm > tol) { 9418839329SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "MatMatSolve: Norm of error %g\n", (double)norm)); 9518839329SStefano Zampini PetscCall(MatConvert(Y, MATAIJ, MAT_INPLACE_MATRIX, &Y)); 9618839329SStefano Zampini PetscCall(MatFilter(Y, PETSC_SMALL, PETSC_TRUE, PETSC_FALSE)); 9718839329SStefano Zampini PetscCall(MatViewFromOptions(Y, NULL, "-aij_err_view")); 9818839329SStefano Zampini } 9918839329SStefano Zampini PetscCall(MatDestroy(&A)); 10018839329SStefano Zampini PetscCall(MatDestroy(&X)); 10118839329SStefano Zampini PetscCall(MatDestroy(&F)); 10218839329SStefano Zampini PetscCall(MatDestroy(&B)); 10318839329SStefano Zampini PetscCall(MatDestroy(&Y)); 10418839329SStefano Zampini PetscCall(PetscFinalize()); 10518839329SStefano Zampini return 0; 10618839329SStefano Zampini } 10718839329SStefano Zampini 10818839329SStefano Zampini /*TEST 10918839329SStefano Zampini 11018839329SStefano Zampini test: 11118839329SStefano Zampini output_file: output/ex62_1.out 11218839329SStefano Zampini suffix: mumps_1 113*7bf14937SSatish Balay requires: datafilespath mumps double !complex !defined(PETSC_USE_64BIT_INDICES) 11418839329SStefano Zampini args: -A ${DATAFILESPATH}/matrices/factorSchur/A.dat -B ${DATAFILESPATH}/matrices/factorSchur/B1.dat -ns {{0 1}} 11518839329SStefano Zampini 11618839329SStefano Zampini test: 11718839329SStefano Zampini output_file: output/ex62_1.out 11818839329SStefano Zampini suffix: mumps_2 119*7bf14937SSatish Balay requires: datafilespath mumps double !complex !defined(PETSC_USE_64BIT_INDICES) 12018839329SStefano Zampini args: -A ${DATAFILESPATH}/matrices/factorSchur/A.dat -B ${DATAFILESPATH}/matrices/factorSchur/B2.dat -ns {{0 1}} 12118839329SStefano Zampini 12218839329SStefano Zampini TEST*/ 123