15e5bbd0aSStefano Zampini #include <petsc/private/pcbddcimpl.h> 25e5bbd0aSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h> 325d06dbeSstefano_zampini #include <petscblaslapack.h> 425d06dbeSstefano_zampini 5d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y) 6d71ae5a4SJacob Faibussowitsch { 725d06dbeSstefano_zampini BDdelta_DN ctx; 825d06dbeSstefano_zampini 925d06dbeSstefano_zampini PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &ctx)); 119566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ctx->BD, x, ctx->work)); 129566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ctx->kBD, ctx->work, y)); 13c0decd05SBarry Smith /* No PC so cannot propagate up failure in KSPSolveTranspose() */ 143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1525d06dbeSstefano_zampini } 1625d06dbeSstefano_zampini 17d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y) 18d71ae5a4SJacob Faibussowitsch { 1925d06dbeSstefano_zampini BDdelta_DN ctx; 2025d06dbeSstefano_zampini 2125d06dbeSstefano_zampini PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &ctx)); 239566063dSJacob Faibussowitsch PetscCall(KSPSolve(ctx->kBD, x, ctx->work)); 24c0decd05SBarry Smith /* No PC so cannot propagate up failure in KSPSolve() */ 259566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->BD, ctx->work, y)); 263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2725d06dbeSstefano_zampini } 2825d06dbeSstefano_zampini 29d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A) 30d71ae5a4SJacob Faibussowitsch { 3125d06dbeSstefano_zampini BDdelta_DN ctx; 3225d06dbeSstefano_zampini 3325d06dbeSstefano_zampini PetscFunctionBegin; 349566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &ctx)); 359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->BD)); 369566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ctx->kBD)); 379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->work)); 389566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4025d06dbeSstefano_zampini } 4125d06dbeSstefano_zampini 42d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx) 43d71ae5a4SJacob Faibussowitsch { 44674ae819SStefano Zampini FETIDPMat_ctx newctx; 45674ae819SStefano Zampini 46674ae819SStefano Zampini PetscFunctionBegin; 479566063dSJacob Faibussowitsch PetscCall(PetscNew(&newctx)); 48674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc)); 50674ae819SStefano Zampini newctx->pc = pc; 51674ae819SStefano Zampini *fetidpmat_ctx = newctx; 523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53674ae819SStefano Zampini } 54674ae819SStefano Zampini 55d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx) 56d71ae5a4SJacob Faibussowitsch { 57674ae819SStefano Zampini FETIDPPC_ctx newctx; 58674ae819SStefano Zampini 59674ae819SStefano Zampini PetscFunctionBegin; 609566063dSJacob Faibussowitsch PetscCall(PetscNew(&newctx)); 61674ae819SStefano Zampini /* increase the reference count for BDDC preconditioner */ 629566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc)); 63674ae819SStefano Zampini newctx->pc = pc; 64674ae819SStefano Zampini *fetidppc_ctx = newctx; 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66674ae819SStefano Zampini } 67674ae819SStefano Zampini 68d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A) 69d71ae5a4SJacob Faibussowitsch { 70674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 71674ae819SStefano Zampini 72674ae819SStefano Zampini PetscFunctionBegin; 739566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &mat_ctx)); 749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->lambda_local)); 759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->temp_solution_D)); 769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->temp_solution_B)); 779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->B_delta)); 789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->B_Ddelta)); 7926699680SStefano Zampini PetscCall(ISDestroy(&mat_ctx->lP_I)); 8026699680SStefano Zampini PetscCall(ISDestroy(&mat_ctx->lP_B)); 819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->B_BB)); 829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->B_BI)); 839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->Bt_BB)); 849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->Bt_BI)); 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_ctx->C)); 869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->rhs_flip)); 879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->vP)); 889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->xPg)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mat_ctx->yPg)); 909566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda)); 919566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda_only)); 929566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&mat_ctx->l2g_p)); 939566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&mat_ctx->g2g_p)); 949566063dSJacob Faibussowitsch PetscCall(PCDestroy(&mat_ctx->pc)); /* decrease PCBDDC reference count */ 959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&mat_ctx->pressure)); 969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&mat_ctx->lagrange)); 979566063dSJacob Faibussowitsch PetscCall(PetscFree(mat_ctx)); 983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 99674ae819SStefano Zampini } 100674ae819SStefano Zampini 101d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) 102d71ae5a4SJacob Faibussowitsch { 103674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 104674ae819SStefano Zampini 105674ae819SStefano Zampini PetscFunctionBegin; 1069566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &pc_ctx)); 1079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc_ctx->lambda_local)); 1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc_ctx->B_Ddelta)); 1099566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&pc_ctx->l2g_lambda)); 1109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc_ctx->S_j)); 1119566063dSJacob Faibussowitsch PetscCall(PCDestroy(&pc_ctx->pc)); /* decrease PCBDDC reference count */ 1129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc_ctx->xPg)); 1139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc_ctx->yPg)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(pc_ctx)); 1153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116674ae819SStefano Zampini } 117674ae819SStefano Zampini 118d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx) 119d71ae5a4SJacob Faibussowitsch { 120674ae819SStefano Zampini PC_IS *pcis = (PC_IS *)fetidpmat_ctx->pc->data; 121674ae819SStefano Zampini PC_BDDC *pcbddc = (PC_BDDC *)fetidpmat_ctx->pc->data; 122674ae819SStefano Zampini PCBDDCGraph mat_graph = pcbddc->mat_graph; 123674ae819SStefano Zampini Mat_IS *matis = (Mat_IS *)fetidpmat_ctx->pc->pmat->data; 124674ae819SStefano Zampini MPI_Comm comm; 12525d06dbeSstefano_zampini Mat ScalingMat, BD1, BD2; 126457d4a33Sstefano_zampini Vec fetidp_global; 127674ae819SStefano Zampini IS IS_l2g_lambda; 128a1c0d0daSstefano_zampini IS subset, subset_mult, subset_n, isvert; 129674ae819SStefano Zampini PetscBool skip_node, fully_redundant; 130674ae819SStefano Zampini PetscInt i, j, k, s, n_boundary_dofs, n_global_lambda, n_vertices, partial_sum; 1316497c311SBarry Smith PetscInt cum, n_local_lambda, n_lambda_for_dof, dual_size, n_neg_values, n_pos_values, buf_size; 13260b1fa21SPierre Jolivet PetscMPIInt rank, size, neigh; 133674ae819SStefano Zampini PetscScalar scalar_value; 134a1c0d0daSstefano_zampini const PetscInt *vertex_indices; 135dc456d91SStefano Zampini PetscInt *dual_dofs_boundary_indices, *aux_local_numbering_1; 136dc456d91SStefano Zampini const PetscInt *aux_global_numbering; 137674ae819SStefano Zampini PetscInt *aux_sums, *cols_B_delta, *l2g_indices; 138674ae819SStefano Zampini PetscScalar *array, *scaling_factors, *vals_B_delta; 13925d06dbeSstefano_zampini PetscScalar **all_factors; 140674ae819SStefano Zampini PetscInt *aux_local_numbering_2; 1419de2952eSStefano Zampini PetscInt *count, **neighbours_set; 142457d4a33Sstefano_zampini PetscLayout llay; 143a1c0d0daSstefano_zampini 144457d4a33Sstefano_zampini /* saddlepoint */ 145457d4a33Sstefano_zampini ISLocalToGlobalMapping l2gmap_p; 146457d4a33Sstefano_zampini PetscLayout play; 147457d4a33Sstefano_zampini IS gP, pP; 148457d4a33Sstefano_zampini PetscInt nPl, nPg, nPgl; 149674ae819SStefano Zampini 150674ae819SStefano Zampini PetscFunctionBegin; 151f4f49eeaSPierre Jolivet PetscCall(PetscObjectGetComm((PetscObject)fetidpmat_ctx->pc, &comm)); 1529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 154674ae819SStefano Zampini 155457d4a33Sstefano_zampini /* saddlepoint */ 156457d4a33Sstefano_zampini nPl = 0; 157457d4a33Sstefano_zampini nPg = 0; 158457d4a33Sstefano_zampini nPgl = 0; 159457d4a33Sstefano_zampini gP = NULL; 160457d4a33Sstefano_zampini pP = NULL; 161457d4a33Sstefano_zampini l2gmap_p = NULL; 162457d4a33Sstefano_zampini play = NULL; 1639566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_pP", (PetscObject *)&pP)); 164022d8d2bSstefano_zampini if (pP) { /* saddle point */ 165457d4a33Sstefano_zampini /* subdomain pressures in global numbering */ 1669566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_gP", (PetscObject *)&gP)); 16728b400f6SJacob Faibussowitsch PetscCheck(gP, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gP not present"); 1689566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(gP, &nPl)); 1699566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->vP)); 1709566063dSJacob Faibussowitsch PetscCall(VecSetSizes(fetidpmat_ctx->vP, nPl, nPl)); 1719566063dSJacob Faibussowitsch PetscCall(VecSetType(fetidpmat_ctx->vP, VECSTANDARD)); 1729566063dSJacob Faibussowitsch PetscCall(VecSetUp(fetidpmat_ctx->vP)); 173457d4a33Sstefano_zampini 174e1214c54Sstefano_zampini /* pressure matrix */ 1759566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_C", (PetscObject *)&fetidpmat_ctx->C)); 176e1214c54Sstefano_zampini if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */ 177457d4a33Sstefano_zampini IS Pg; 178457d4a33Sstefano_zampini 1799566063dSJacob Faibussowitsch PetscCall(ISRenumber(gP, NULL, &nPg, &Pg)); 1809566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(Pg, &l2gmap_p)); 1819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&Pg)); 1829566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &play)); 1839566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(play, 1)); 1849566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(play, nPg)); 1859566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pP, &nPgl)); 1869566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(play, nPgl)); 1879566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(play)); 188457d4a33Sstefano_zampini } else { 1899566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->C)); 1909566063dSJacob Faibussowitsch PetscCall(MatISGetLocalToGlobalMapping(fetidpmat_ctx->C, &l2gmap_p, NULL)); 1919566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)l2gmap_p)); 1929566063dSJacob Faibussowitsch PetscCall(MatGetSize(fetidpmat_ctx->C, &nPg, NULL)); 1939566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(fetidpmat_ctx->C, NULL, &nPgl)); 1949566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(fetidpmat_ctx->C, NULL, &llay)); 1959566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(llay, &play)); 196457d4a33Sstefano_zampini } 1979566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->xPg)); 1989566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->yPg)); 199457d4a33Sstefano_zampini 200e1214c54Sstefano_zampini /* import matrices for pressures coupling */ 2019566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BI", (PetscObject *)&fetidpmat_ctx->B_BI)); 20228b400f6SJacob Faibussowitsch PetscCheck(fetidpmat_ctx->B_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BI not present"); 2039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI)); 204457d4a33Sstefano_zampini 2059566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BB", (PetscObject *)&fetidpmat_ctx->B_BB)); 20628b400f6SJacob Faibussowitsch PetscCheck(fetidpmat_ctx->B_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BB not present"); 2079566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB)); 208457d4a33Sstefano_zampini 2099566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BI", (PetscObject *)&fetidpmat_ctx->Bt_BI)); 21028b400f6SJacob Faibussowitsch PetscCheck(fetidpmat_ctx->Bt_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BI not present"); 2119566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI)); 212457d4a33Sstefano_zampini 2139566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BB", (PetscObject *)&fetidpmat_ctx->Bt_BB)); 21428b400f6SJacob Faibussowitsch PetscCheck(fetidpmat_ctx->Bt_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BB not present"); 2159566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB)); 2166cc1294bSstefano_zampini 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_flip", (PetscObject *)&fetidpmat_ctx->rhs_flip)); 2181baa6e33SBarry Smith if (fetidpmat_ctx->rhs_flip) PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip)); 21926699680SStefano Zampini 22026699680SStefano Zampini PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_lP_I", (PetscObject *)&fetidpmat_ctx->lP_I)); 22126699680SStefano Zampini PetscCheck(fetidpmat_ctx->lP_I, PETSC_COMM_SELF, PETSC_ERR_PLIB, "lP_I not present"); 22226699680SStefano Zampini PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->lP_I)); 22326699680SStefano Zampini 22426699680SStefano Zampini PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_lP_B", (PetscObject *)&fetidpmat_ctx->lP_B)); 22526699680SStefano Zampini PetscCheck(fetidpmat_ctx->lP_B, PETSC_COMM_SELF, PETSC_ERR_PLIB, "lP_B not present"); 22626699680SStefano Zampini PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->lP_B)); 227457d4a33Sstefano_zampini } 228457d4a33Sstefano_zampini 229674ae819SStefano Zampini /* Default type of lagrange multipliers is non-redundant */ 230329cd39dSStefano Zampini fully_redundant = fetidpmat_ctx->fully_redundant; 231674ae819SStefano Zampini 232674ae819SStefano Zampini /* Evaluate local and global number of lagrange multipliers */ 2339566063dSJacob Faibussowitsch PetscCall(VecSet(pcis->vec1_N, 0.0)); 234674ae819SStefano Zampini n_local_lambda = 0; 235674ae819SStefano Zampini partial_sum = 0; 236674ae819SStefano Zampini n_boundary_dofs = 0; 237a1c0d0daSstefano_zampini 238674ae819SStefano Zampini /* Get Vertices used to define the BDDC */ 2399de2952eSStefano Zampini PetscCall(PCBDDCGraphGetCandidatesIS(mat_graph, NULL, NULL, NULL, NULL, &isvert)); 2409566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isvert, &n_vertices)); 2419566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isvert, &vertex_indices)); 242a1c0d0daSstefano_zampini 243674ae819SStefano Zampini dual_size = pcis->n_B - n_vertices; 2449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dual_size, &dual_dofs_boundary_indices)); 2459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_1)); 2469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_2)); 247674ae819SStefano Zampini 2489de2952eSStefano Zampini /* the code below does not support multiple subdomains per process 2499de2952eSStefano Zampini error out in this case 2509de2952eSStefano Zampini TODO: I guess I can use PetscSFGetMultiSF and the code will be easier and more general */ 2519de2952eSStefano Zampini PetscCall(PetscMalloc2(pcis->n, &count, pcis->n, &neighbours_set)); 2529de2952eSStefano Zampini for (i = 0, j = 0; i < pcis->n; i++) j += mat_graph->nodes[i].count; 2539de2952eSStefano Zampini if (pcis->n) PetscCall(PetscMalloc1(j, &neighbours_set[0])); 254674ae819SStefano Zampini for (i = 0; i < pcis->n; i++) { 2559de2952eSStefano Zampini PCBDDCGraphNode *node = &mat_graph->nodes[i]; 2569de2952eSStefano Zampini 2579de2952eSStefano Zampini count[i] = 0; 2589de2952eSStefano Zampini for (j = 0; j < node->count; j++) { 2599de2952eSStefano Zampini if (node->neighbours_set[j] == rank) continue; 2609de2952eSStefano Zampini neighbours_set[i][count[i]++] = node->neighbours_set[j]; 2619de2952eSStefano Zampini } 2629de2952eSStefano Zampini PetscCheck(count[i] == node->count - 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Multiple subdomains per process not supported"); 2639de2952eSStefano Zampini s = count[i]; 2649de2952eSStefano Zampini PetscCall(PetscSortRemoveDupsInt(count + i, neighbours_set[i])); 2659de2952eSStefano Zampini PetscCheck(s == count[i], PETSC_COMM_SELF, PETSC_ERR_SUP, "Multiple subdomains per process not supported"); 2669de2952eSStefano Zampini if (i != pcis->n - 1) neighbours_set[i + 1] = neighbours_set[i] + count[i]; 2679de2952eSStefano Zampini } 2689de2952eSStefano Zampini 2699de2952eSStefano Zampini PetscCall(VecGetArray(pcis->vec1_N, &array)); 2709de2952eSStefano Zampini for (i = 0, s = 0; i < pcis->n; i++) { 2719de2952eSStefano Zampini j = count[i]; /* RECALL: count[i] does not count myself */ 2722d456bbaSstefano_zampini if (j > 0) n_boundary_dofs++; 273674ae819SStefano Zampini skip_node = PETSC_FALSE; 274674ae819SStefano Zampini if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */ 275674ae819SStefano Zampini skip_node = PETSC_TRUE; 276674ae819SStefano Zampini s++; 277674ae819SStefano Zampini } 2782d456bbaSstefano_zampini if (j < 1) skip_node = PETSC_TRUE; 2799de2952eSStefano Zampini if (mat_graph->nodes[i].special_dof == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE; 280674ae819SStefano Zampini if (!skip_node) { 281674ae819SStefano Zampini if (fully_redundant) { 282674ae819SStefano Zampini /* fully redundant set of lagrange multipliers */ 283674ae819SStefano Zampini n_lambda_for_dof = (j * (j + 1)) / 2; 284674ae819SStefano Zampini } else { 285674ae819SStefano Zampini n_lambda_for_dof = j; 286674ae819SStefano Zampini } 287674ae819SStefano Zampini n_local_lambda += j; 288674ae819SStefano Zampini /* needed to evaluate global number of lagrange multipliers */ 289674ae819SStefano Zampini array[i] = (1.0 * n_lambda_for_dof) / (j + 1.0); /* already scaled for the next global sum */ 290674ae819SStefano Zampini /* store some data needed */ 291674ae819SStefano Zampini dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs - 1; 292674ae819SStefano Zampini aux_local_numbering_1[partial_sum] = i; 293674ae819SStefano Zampini aux_local_numbering_2[partial_sum] = n_lambda_for_dof; 294674ae819SStefano Zampini partial_sum++; 295674ae819SStefano Zampini } 296674ae819SStefano Zampini } 2979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(pcis->vec1_N, &array)); 2989566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isvert, &vertex_indices)); 2999de2952eSStefano Zampini PetscCall(PCBDDCGraphRestoreCandidatesIS(mat_graph, NULL, NULL, NULL, NULL, &isvert)); 3002d456bbaSstefano_zampini dual_size = partial_sum; 301674ae819SStefano Zampini 302674ae819SStefano Zampini /* compute global ordering of lagrange multipliers and associate l2g map */ 3039566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_1, PETSC_COPY_VALUES, &subset_n)); 3049566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingApplyIS(pcis->mapping, subset_n, &subset)); 3059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&subset_n)); 3069566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_2, PETSC_OWN_POINTER, &subset_mult)); 3079566063dSJacob Faibussowitsch PetscCall(ISRenumber(subset, subset_mult, &fetidpmat_ctx->n_lambda, &subset_n)); 3089566063dSJacob Faibussowitsch PetscCall(ISDestroy(&subset)); 3093d996552SStefano Zampini 31076bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 3119566063dSJacob Faibussowitsch PetscCall(VecSet(pcis->vec1_global, 0.0)); 3129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE)); 3139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE)); 3149566063dSJacob Faibussowitsch PetscCall(VecSum(pcis->vec1_global, &scalar_value)); 3154fe826edSStefano Zampini i = (PetscInt)PetscRealPart(scalar_value); 31663a3b9bcSJacob Faibussowitsch PetscCheck(i == fetidpmat_ctx->n_lambda, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Global number of multipliers mismatch! (%" PetscInt_FMT " != %" PetscInt_FMT ")", fetidpmat_ctx->n_lambda, i); 31776bd3646SJed Brown } 318674ae819SStefano Zampini 319674ae819SStefano Zampini /* init data for scaling factors exchange */ 32025d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 32125d06dbeSstefano_zampini PetscInt *ptrs_buffer, neigh_position; 32225d06dbeSstefano_zampini PetscScalar *send_buffer, *recv_buffer; 32325d06dbeSstefano_zampini MPI_Request *send_reqs, *recv_reqs; 324835f2295SStefano Zampini PetscMPIInt nreqs; 32525d06dbeSstefano_zampini 326674ae819SStefano Zampini partial_sum = 0; 3279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pcis->n_neigh, &ptrs_buffer)); 3289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &send_reqs)); 3299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &recv_reqs)); 3309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pcis->n + 1, &all_factors)); 3314b2aedd3SStefano Zampini if (pcis->n_neigh > 0) ptrs_buffer[0] = 0; 332674ae819SStefano Zampini for (i = 1; i < pcis->n_neigh; i++) { 333674ae819SStefano Zampini partial_sum += pcis->n_shared[i]; 334674ae819SStefano Zampini ptrs_buffer[i] = ptrs_buffer[i - 1] + pcis->n_shared[i]; 335674ae819SStefano Zampini } 3369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(partial_sum, &send_buffer)); 3379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(partial_sum, &recv_buffer)); 3389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(partial_sum, &all_factors[0])); 339674ae819SStefano Zampini for (i = 0; i < pcis->n - 1; i++) { 3409de2952eSStefano Zampini j = count[i]; 341674ae819SStefano Zampini all_factors[i + 1] = all_factors[i] + j; 342674ae819SStefano Zampini } 34325d06dbeSstefano_zampini 344674ae819SStefano Zampini /* scatter B scaling to N vec */ 3459566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 3469566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 347674ae819SStefano Zampini /* communications */ 3489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pcis->vec1_N, (const PetscScalar **)&array)); 349674ae819SStefano Zampini for (i = 1; i < pcis->n_neigh; i++) { 350ad540459SPierre Jolivet for (j = 0; j < pcis->n_shared[i]; j++) send_buffer[ptrs_buffer[i - 1] + j] = array[pcis->shared[i][j]]; 3516497c311SBarry Smith buf_size = ptrs_buffer[i] - ptrs_buffer[i - 1]; 35260b1fa21SPierre Jolivet PetscCall(PetscMPIIntCast(pcis->neigh[i], &neigh)); 35360b1fa21SPierre Jolivet PetscCallMPI(MPIU_Isend(&send_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &send_reqs[i - 1])); 35460b1fa21SPierre Jolivet PetscCallMPI(MPIU_Irecv(&recv_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &recv_reqs[i - 1])); 355674ae819SStefano Zampini } 3569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pcis->vec1_N, (const PetscScalar **)&array)); 357835f2295SStefano Zampini PetscCall(PetscMPIIntCast(pcis->n_neigh - 1, &nreqs)); 358835f2295SStefano Zampini if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(nreqs, recv_reqs, MPI_STATUSES_IGNORE)); 359674ae819SStefano Zampini /* put values in correct places */ 360674ae819SStefano Zampini for (i = 1; i < pcis->n_neigh; i++) { 361674ae819SStefano Zampini for (j = 0; j < pcis->n_shared[i]; j++) { 362674ae819SStefano Zampini k = pcis->shared[i][j]; 363674ae819SStefano Zampini neigh_position = 0; 364ac530a7eSPierre Jolivet while (neighbours_set[k][neigh_position] != pcis->neigh[i]) neigh_position++; 365674ae819SStefano Zampini all_factors[k][neigh_position] = recv_buffer[ptrs_buffer[i - 1] + j]; 366674ae819SStefano Zampini } 367674ae819SStefano Zampini } 368835f2295SStefano Zampini if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(nreqs, send_reqs, MPI_STATUSES_IGNORE)); 3699566063dSJacob Faibussowitsch PetscCall(PetscFree(send_reqs)); 3709566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_reqs)); 3719566063dSJacob Faibussowitsch PetscCall(PetscFree(send_buffer)); 3729566063dSJacob Faibussowitsch PetscCall(PetscFree(recv_buffer)); 3739566063dSJacob Faibussowitsch PetscCall(PetscFree(ptrs_buffer)); 37425d06dbeSstefano_zampini } 375674ae819SStefano Zampini 376674ae819SStefano Zampini /* Compute B and B_delta (local actions) */ 3779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pcis->n_neigh, &aux_sums)); 3789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local_lambda, &l2g_indices)); 3799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local_lambda, &vals_B_delta)); 3809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local_lambda, &cols_B_delta)); 38125d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 3829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local_lambda, &scaling_factors)); 38325d06dbeSstefano_zampini } else { 38425d06dbeSstefano_zampini scaling_factors = NULL; 38525d06dbeSstefano_zampini all_factors = NULL; 38625d06dbeSstefano_zampini } 3879566063dSJacob Faibussowitsch PetscCall(ISGetIndices(subset_n, &aux_global_numbering)); 388674ae819SStefano Zampini partial_sum = 0; 389dc456d91SStefano Zampini cum = 0; 390674ae819SStefano Zampini for (i = 0; i < dual_size; i++) { 391dc456d91SStefano Zampini n_global_lambda = aux_global_numbering[cum]; 3929de2952eSStefano Zampini j = count[aux_local_numbering_1[i]]; 393674ae819SStefano Zampini aux_sums[0] = 0; 394ad540459SPierre Jolivet for (s = 1; s < j; s++) aux_sums[s] = aux_sums[s - 1] + j - s + 1; 39525d06dbeSstefano_zampini if (all_factors) array = all_factors[aux_local_numbering_1[i]]; 396674ae819SStefano Zampini n_neg_values = 0; 397ac530a7eSPierre Jolivet while (n_neg_values < j && neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) n_neg_values++; 398674ae819SStefano Zampini n_pos_values = j - n_neg_values; 399674ae819SStefano Zampini if (fully_redundant) { 400674ae819SStefano Zampini for (s = 0; s < n_neg_values; s++) { 401674ae819SStefano Zampini l2g_indices[partial_sum + s] = aux_sums[s] + n_neg_values - s - 1 + n_global_lambda; 402674ae819SStefano Zampini cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i]; 403674ae819SStefano Zampini vals_B_delta[partial_sum + s] = -1.0; 40425d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s] = array[s]; 405674ae819SStefano Zampini } 406674ae819SStefano Zampini for (s = 0; s < n_pos_values; s++) { 407674ae819SStefano Zampini l2g_indices[partial_sum + s + n_neg_values] = aux_sums[n_neg_values] + s + n_global_lambda; 408674ae819SStefano Zampini cols_B_delta[partial_sum + s + n_neg_values] = dual_dofs_boundary_indices[i]; 409674ae819SStefano Zampini vals_B_delta[partial_sum + s + n_neg_values] = 1.0; 41025d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s + n_neg_values] = array[s + n_neg_values]; 411674ae819SStefano Zampini } 412674ae819SStefano Zampini partial_sum += j; 413674ae819SStefano Zampini } else { 414674ae819SStefano Zampini /* l2g_indices and default cols and vals of B_delta */ 415674ae819SStefano Zampini for (s = 0; s < j; s++) { 416674ae819SStefano Zampini l2g_indices[partial_sum + s] = n_global_lambda + s; 417674ae819SStefano Zampini cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i]; 418674ae819SStefano Zampini vals_B_delta[partial_sum + s] = 0.0; 419674ae819SStefano Zampini } 420674ae819SStefano Zampini /* B_delta */ 421674ae819SStefano Zampini if (n_neg_values > 0) { /* there's a rank next to me to the left */ 422674ae819SStefano Zampini vals_B_delta[partial_sum + n_neg_values - 1] = -1.0; 423674ae819SStefano Zampini } 424674ae819SStefano Zampini if (n_neg_values < j) { /* there's a rank next to me to the right */ 425674ae819SStefano Zampini vals_B_delta[partial_sum + n_neg_values] = 1.0; 426674ae819SStefano Zampini } 427674ae819SStefano Zampini /* scaling as in Klawonn-Widlund 1999 */ 42825d06dbeSstefano_zampini if (!pcbddc->use_deluxe_scaling) { 429674ae819SStefano Zampini for (s = 0; s < n_neg_values; s++) { 430674ae819SStefano Zampini scalar_value = 0.0; 43125d06dbeSstefano_zampini for (k = 0; k < s + 1; k++) scalar_value += array[k]; 432674ae819SStefano Zampini scaling_factors[partial_sum + s] = -scalar_value; 433674ae819SStefano Zampini } 434674ae819SStefano Zampini for (s = 0; s < n_pos_values; s++) { 435674ae819SStefano Zampini scalar_value = 0.0; 43625d06dbeSstefano_zampini for (k = s + n_neg_values; k < j; k++) scalar_value += array[k]; 437674ae819SStefano Zampini scaling_factors[partial_sum + s + n_neg_values] = scalar_value; 438674ae819SStefano Zampini } 43925d06dbeSstefano_zampini } 440674ae819SStefano Zampini partial_sum += j; 441674ae819SStefano Zampini } 442dc456d91SStefano Zampini cum += aux_local_numbering_2[i]; 443674ae819SStefano Zampini } 4449566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(subset_n, &aux_global_numbering)); 4459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&subset_mult)); 4469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&subset_n)); 4479566063dSJacob Faibussowitsch PetscCall(PetscFree(aux_sums)); 4489566063dSJacob Faibussowitsch PetscCall(PetscFree(aux_local_numbering_1)); 4499566063dSJacob Faibussowitsch PetscCall(PetscFree(dual_dofs_boundary_indices)); 45025d06dbeSstefano_zampini if (all_factors) { 4519566063dSJacob Faibussowitsch PetscCall(PetscFree(all_factors[0])); 4529566063dSJacob Faibussowitsch PetscCall(PetscFree(all_factors)); 45325d06dbeSstefano_zampini } 4549de2952eSStefano Zampini if (pcis->n) PetscCall(PetscFree(neighbours_set[0])); 4559de2952eSStefano Zampini PetscCall(PetscFree2(count, neighbours_set)); 456674ae819SStefano Zampini 457674ae819SStefano Zampini /* Create local part of B_delta */ 4589566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_delta)); 4599566063dSJacob Faibussowitsch PetscCall(MatSetSizes(fetidpmat_ctx->B_delta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B)); 4609566063dSJacob Faibussowitsch PetscCall(MatSetType(fetidpmat_ctx->B_delta, MATSEQAIJ)); 4619566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta, 1, NULL)); 4629566063dSJacob Faibussowitsch PetscCall(MatSetOption(fetidpmat_ctx->B_delta, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 46348a46eb9SPierre Jolivet for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(fetidpmat_ctx->B_delta, i, cols_B_delta[i], vals_B_delta[i], INSERT_VALUES)); 4649566063dSJacob Faibussowitsch PetscCall(PetscFree(vals_B_delta)); 4659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY)); 4669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY)); 467674ae819SStefano Zampini 46825d06dbeSstefano_zampini BD1 = NULL; 46925d06dbeSstefano_zampini BD2 = NULL; 470674ae819SStefano Zampini if (fully_redundant) { 47128b400f6SJacob Faibussowitsch PetscCheck(!pcbddc->use_deluxe_scaling, comm, PETSC_ERR_SUP, "Deluxe FETIDP with fully-redundant multipliers to be implemented"); 4729566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &ScalingMat)); 4739566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ScalingMat, n_local_lambda, n_local_lambda, n_local_lambda, n_local_lambda)); 4749566063dSJacob Faibussowitsch PetscCall(MatSetType(ScalingMat, MATSEQAIJ)); 4759566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ScalingMat, 1, NULL)); 47648a46eb9SPierre Jolivet for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(ScalingMat, i, i, scaling_factors[i], INSERT_VALUES)); 4779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ScalingMat, MAT_FINAL_ASSEMBLY)); 4789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ScalingMat, MAT_FINAL_ASSEMBLY)); 479fb842aefSJose E. Roman PetscCall(MatMatMult(ScalingMat, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &fetidpmat_ctx->B_Ddelta)); 4809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ScalingMat)); 481674ae819SStefano Zampini } else { 4829566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_Ddelta)); 4839566063dSJacob Faibussowitsch PetscCall(MatSetSizes(fetidpmat_ctx->B_Ddelta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B)); 484524221abSStefano Zampini if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) { 4859566063dSJacob Faibussowitsch PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSEQAIJ)); 4869566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta, 1, NULL)); 48748a46eb9SPierre Jolivet for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(fetidpmat_ctx->B_Ddelta, i, cols_B_delta[i], scaling_factors[i], INSERT_VALUES)); 4889566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY)); 4899566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY)); 49025d06dbeSstefano_zampini } else { 49125d06dbeSstefano_zampini /* scaling as in Klawonn-Widlund 1999 */ 49225d06dbeSstefano_zampini PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx; 49325d06dbeSstefano_zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 49425d06dbeSstefano_zampini Mat T; 49525d06dbeSstefano_zampini PetscScalar *W, lwork, *Bwork; 496451a39c7SStefano Zampini const PetscInt *idxs = NULL; 49725d06dbeSstefano_zampini PetscInt cum, mss, *nnz; 49825d06dbeSstefano_zampini PetscBLASInt *pivots, B_lwork, B_N, B_ierr; 49925d06dbeSstefano_zampini 50028b400f6SJacob Faibussowitsch PetscCheck(pcbddc->deluxe_singlemat, comm, PETSC_ERR_USER, "Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat"); 50125d06dbeSstefano_zampini mss = 0; 5029566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(pcis->n_B, &nnz)); 50325d06dbeSstefano_zampini if (sub_schurs->is_Ej_all) { 5049566063dSJacob Faibussowitsch PetscCall(ISGetIndices(sub_schurs->is_Ej_all, &idxs)); 50525d06dbeSstefano_zampini for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) { 50625d06dbeSstefano_zampini PetscInt subset_size; 50725d06dbeSstefano_zampini 5089566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 50925d06dbeSstefano_zampini for (j = 0; j < subset_size; j++) nnz[idxs[j + cum]] = subset_size; 51025d06dbeSstefano_zampini mss = PetscMax(mss, subset_size); 51125d06dbeSstefano_zampini cum += subset_size; 51225d06dbeSstefano_zampini } 51325d06dbeSstefano_zampini } 5149566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &T)); 5159566063dSJacob Faibussowitsch PetscCall(MatSetSizes(T, pcis->n_B, pcis->n_B, pcis->n_B, pcis->n_B)); 5169566063dSJacob Faibussowitsch PetscCall(MatSetType(T, MATSEQAIJ)); 5179566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(T, 0, nnz)); 5189566063dSJacob Faibussowitsch PetscCall(PetscFree(nnz)); 51925d06dbeSstefano_zampini 52025d06dbeSstefano_zampini /* workspace allocation */ 5217f00194aSstefano_zampini B_lwork = 0; 5227f00194aSstefano_zampini if (mss) { 523be8bcbd0Sstefano_zampini PetscScalar dummy = 1; 524be8bcbd0Sstefano_zampini 52525d06dbeSstefano_zampini B_lwork = -1; 5269566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(mss, &B_N)); 5279566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 528792fecdfSBarry Smith PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, &dummy, &B_N, &B_N, &lwork, &B_lwork, &B_ierr)); 5299566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 530835f2295SStefano Zampini PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to GETRI Lapack routine %" PetscBLASInt_FMT, B_ierr); 5319566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork), &B_lwork)); 5327f00194aSstefano_zampini } 5339566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(mss * mss, &W, mss, &pivots, B_lwork, &Bwork)); 53425d06dbeSstefano_zampini 53525d06dbeSstefano_zampini for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) { 5361683a169SBarry Smith const PetscScalar *M; 53725d06dbeSstefano_zampini PetscInt subset_size; 53825d06dbeSstefano_zampini 5399566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 5409566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(subset_size, &B_N)); 5419566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(deluxe_ctx->seq_mat[i], &M)); 5429566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(W, M, subset_size * subset_size)); 5439566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i], &M)); 5449566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 545792fecdfSBarry Smith PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, W, &B_N, pivots, &B_ierr)); 546835f2295SStefano Zampini PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %" PetscBLASInt_FMT, B_ierr); 547792fecdfSBarry Smith PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, W, &B_N, pivots, Bwork, &B_lwork, &B_ierr)); 548835f2295SStefano Zampini PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %" PetscBLASInt_FMT, B_ierr); 5499566063dSJacob Faibussowitsch PetscCall(PetscFPTrapPop()); 550451a39c7SStefano Zampini /* silent static analyzer */ 55128b400f6SJacob Faibussowitsch PetscCheck(idxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IDXS not present"); 5529566063dSJacob Faibussowitsch PetscCall(MatSetValues(T, subset_size, idxs + cum, subset_size, idxs + cum, W, INSERT_VALUES)); 55325d06dbeSstefano_zampini cum += subset_size; 55425d06dbeSstefano_zampini } 5559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY)); 5569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY)); 557fb842aefSJose E. Roman PetscCall(MatMatTransposeMult(T, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &BD1)); 558fb842aefSJose E. Roman PetscCall(MatMatMult(fetidpmat_ctx->B_delta, BD1, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &BD2)); 5599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 5609566063dSJacob Faibussowitsch PetscCall(PetscFree3(W, pivots, Bwork)); 56148a46eb9SPierre Jolivet if (sub_schurs->is_Ej_all) PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all, &idxs)); 56225d06dbeSstefano_zampini } 563674ae819SStefano Zampini } 5649566063dSJacob Faibussowitsch PetscCall(PetscFree(scaling_factors)); 5659566063dSJacob Faibussowitsch PetscCall(PetscFree(cols_B_delta)); 566674ae819SStefano Zampini 567457d4a33Sstefano_zampini /* Layout of multipliers */ 5689566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &llay)); 5699566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(llay, 1)); 5709566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetSize(llay, fetidpmat_ctx->n_lambda)); 5719566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(llay)); 5729566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(llay, &fetidpmat_ctx->n)); 573457d4a33Sstefano_zampini 574457d4a33Sstefano_zampini /* Local work vector of multipliers */ 5759566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->lambda_local)); 5769566063dSJacob Faibussowitsch PetscCall(VecSetSizes(fetidpmat_ctx->lambda_local, n_local_lambda, n_local_lambda)); 5779566063dSJacob Faibussowitsch PetscCall(VecSetType(fetidpmat_ctx->lambda_local, VECSEQ)); 578457d4a33Sstefano_zampini 57925d06dbeSstefano_zampini if (BD2) { 58025d06dbeSstefano_zampini ISLocalToGlobalMapping l2g; 58125d06dbeSstefano_zampini Mat T, TA, *pT; 58225d06dbeSstefano_zampini IS is; 58325d06dbeSstefano_zampini PetscInt nl, N; 58425d06dbeSstefano_zampini BDdelta_DN ctx; 58525d06dbeSstefano_zampini 5869566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetLocalSize(llay, &nl)); 5879566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(llay, &N)); 5889566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &T)); 5899566063dSJacob Faibussowitsch PetscCall(MatSetSizes(T, nl, nl, N, N)); 5909566063dSJacob Faibussowitsch PetscCall(MatSetType(T, MATIS)); 5919566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &l2g)); 5929566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(T, l2g, l2g)); 5939566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&l2g)); 5949566063dSJacob Faibussowitsch PetscCall(MatISSetLocalMat(T, BD2)); 5959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY)); 5969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY)); 5979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BD2)); 5989566063dSJacob Faibussowitsch PetscCall(MatConvert(T, MATAIJ, MAT_INITIAL_MATRIX, &TA)); 5999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 6009566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_USE_POINTER, &is)); 6019566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(TA, 1, &is, &is, MAT_INITIAL_MATRIX, &pT)); 6029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&TA)); 6039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 60425d06dbeSstefano_zampini BD2 = pT[0]; 6059566063dSJacob Faibussowitsch PetscCall(PetscFree(pT)); 60625d06dbeSstefano_zampini 60725d06dbeSstefano_zampini /* B_Ddelta for non-redundant multipliers with deluxe scaling */ 6089566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 6099566063dSJacob Faibussowitsch PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSHELL)); 6109566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(fetidpmat_ctx->B_Ddelta, ctx)); 611*57d50842SBarry Smith PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT, (PetscErrorCodeFn *)MatMult_BDdelta_deluxe_nonred)); 612*57d50842SBarry Smith PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)MatMultTranspose_BDdelta_deluxe_nonred)); 613*57d50842SBarry Smith PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_DESTROY, (PetscErrorCodeFn *)MatDestroy_BDdelta_deluxe_nonred)); 6149566063dSJacob Faibussowitsch PetscCall(MatSetUp(fetidpmat_ctx->B_Ddelta)); 61525d06dbeSstefano_zampini 6169566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)BD1)); 61725d06dbeSstefano_zampini ctx->BD = BD1; 6189566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->kBD)); 6193821be0aSBarry Smith PetscCall(KSPSetNestLevel(ctx->kBD, fetidpmat_ctx->pc->kspnestlevel)); 6209566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ctx->kBD, BD2, BD2)); 6219566063dSJacob Faibussowitsch PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local, &ctx->work)); 62225d06dbeSstefano_zampini fetidpmat_ctx->deluxe_nonred = PETSC_TRUE; 62325d06dbeSstefano_zampini } 6249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BD1)); 6259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&BD2)); 62625d06dbeSstefano_zampini 62725d06dbeSstefano_zampini /* fetidpmat sizes */ 62825d06dbeSstefano_zampini fetidpmat_ctx->n += nPgl; 62925d06dbeSstefano_zampini fetidpmat_ctx->N = fetidpmat_ctx->n_lambda + nPg; 63025d06dbeSstefano_zampini 631457d4a33Sstefano_zampini /* Global vector for FETI-DP linear system */ 6329566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &fetidp_global)); 6339566063dSJacob Faibussowitsch PetscCall(VecSetSizes(fetidp_global, fetidpmat_ctx->n, fetidpmat_ctx->N)); 6349566063dSJacob Faibussowitsch PetscCall(VecSetType(fetidp_global, VECMPI)); 6359566063dSJacob Faibussowitsch PetscCall(VecSetUp(fetidp_global)); 636457d4a33Sstefano_zampini 6379eec4de8Sstefano_zampini /* Decide layout for fetidp dofs: if it is a saddle point problem 6389eec4de8Sstefano_zampini pressure is ordered first in the local part of the global vector 639457d4a33Sstefano_zampini of the FETI-DP linear system */ 640457d4a33Sstefano_zampini if (nPg) { 641e1214c54Sstefano_zampini Vec v; 642af140850Sstefano_zampini IS IS_l2g_p, ais; 643457d4a33Sstefano_zampini PetscLayout alay; 644457d4a33Sstefano_zampini const PetscInt *idxs, *pranges, *aranges, *lranges; 645af140850Sstefano_zampini PetscInt *l2g_indices_p, rst; 646e1214c54Sstefano_zampini PetscMPIInt rank; 647457d4a33Sstefano_zampini 6489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nPl, &l2g_indices_p)); 6499566063dSJacob Faibussowitsch PetscCall(VecGetLayout(fetidp_global, &alay)); 6509566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(alay, &aranges)); 6519566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(play, &pranges)); 6529566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(llay, &lranges)); 653e1214c54Sstefano_zampini 6549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global), &rank)); 6559566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), pranges[rank + 1] - pranges[rank], aranges[rank], 1, &fetidpmat_ctx->pressure)); 6569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure, "F_P")); 6579566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), lranges[rank + 1] - lranges[rank], aranges[rank] + pranges[rank + 1] - pranges[rank], 1, &fetidpmat_ctx->lagrange)); 6589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange, "F_L")); 6599566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(l2gmap_p, &idxs)); 660af140850Sstefano_zampini /* shift local to global indices for pressure */ 661457d4a33Sstefano_zampini for (i = 0; i < nPl; i++) { 662131c27b5Sprj- PetscMPIInt owner; 663457d4a33Sstefano_zampini 6649566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(play, idxs[i], &owner)); 665457d4a33Sstefano_zampini l2g_indices_p[i] = idxs[i] - pranges[owner] + aranges[owner]; 666457d4a33Sstefano_zampini } 6679566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(l2gmap_p, &idxs)); 6689566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nPl, l2g_indices_p, PETSC_OWN_POINTER, &IS_l2g_p)); 669af140850Sstefano_zampini 670e1214c54Sstefano_zampini /* local to global scatter for pressure */ 6719566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(fetidpmat_ctx->vP, NULL, fetidp_global, IS_l2g_p, &fetidpmat_ctx->l2g_p)); 6729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&IS_l2g_p)); 673457d4a33Sstefano_zampini 674e1214c54Sstefano_zampini /* scatter for lagrange multipliers only */ 6759566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &v)); 6769566063dSJacob Faibussowitsch PetscCall(VecSetType(v, VECSTANDARD)); 6779566063dSJacob Faibussowitsch PetscCall(VecSetLayout(v, llay)); 6789566063dSJacob Faibussowitsch PetscCall(VecSetUp(v)); 6799566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &ais)); 6809566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, v, ais, &fetidpmat_ctx->l2g_lambda_only)); 6819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ais)); 6829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 683e1214c54Sstefano_zampini 684af140850Sstefano_zampini /* shift local to global indices for multipliers */ 685457d4a33Sstefano_zampini for (i = 0; i < n_local_lambda; i++) { 686131c27b5Sprj- PetscInt ps; 687131c27b5Sprj- PetscMPIInt owner; 688457d4a33Sstefano_zampini 6899566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(llay, l2g_indices[i], &owner)); 690457d4a33Sstefano_zampini ps = pranges[owner + 1] - pranges[owner]; 691457d4a33Sstefano_zampini l2g_indices[i] = l2g_indices[i] - lranges[owner] + aranges[owner] + ps; 692457d4a33Sstefano_zampini } 693457d4a33Sstefano_zampini 694e1214c54Sstefano_zampini /* scatter from alldofs to pressures global fetidp vector */ 6959566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(alay, &rst, NULL)); 6969566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, nPgl, rst, 1, &ais)); 6979566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(pcis->vec1_global, pP, fetidp_global, ais, &fetidpmat_ctx->g2g_p)); 6989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ais)); 699457d4a33Sstefano_zampini } 7009566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&llay)); 7019566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&play)); 7029566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_OWN_POINTER, &IS_l2g_lambda)); 703a1c0d0daSstefano_zampini 7049cc7774eSstefano_zampini /* scatter from local to global multipliers */ 7059566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, fetidp_global, IS_l2g_lambda, &fetidpmat_ctx->l2g_lambda)); 7069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&IS_l2g_lambda)); 7079566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&l2gmap_p)); 7089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fetidp_global)); 709457d4a33Sstefano_zampini 710a1c0d0daSstefano_zampini /* Create some work vectors needed by fetidp */ 7119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pcis->vec1_B, &fetidpmat_ctx->temp_solution_B)); 7129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(pcis->vec1_D, &fetidpmat_ctx->temp_solution_D)); 7133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 714674ae819SStefano Zampini } 715674ae819SStefano Zampini 716d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) 717d71ae5a4SJacob Faibussowitsch { 718674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 71925d06dbeSstefano_zampini PC_BDDC *pcbddc = (PC_BDDC *)fetidppc_ctx->pc->data; 72025d06dbeSstefano_zampini PC_IS *pcis = (PC_IS *)fetidppc_ctx->pc->data; 721f28b6018SStefano Zampini PetscBool lumped = PETSC_FALSE; 722674ae819SStefano Zampini 723674ae819SStefano Zampini PetscFunctionBegin; 7249566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(fetimat, &mat_ctx)); 725674ae819SStefano Zampini /* get references from objects created when setting up feti mat context */ 7269566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->lambda_local)); 727674ae819SStefano Zampini fetidppc_ctx->lambda_local = mat_ctx->lambda_local; 7289566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->B_Ddelta)); 729674ae819SStefano Zampini fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta; 73025d06dbeSstefano_zampini if (mat_ctx->deluxe_nonred) { 73125d06dbeSstefano_zampini PC pc, mpc; 73225d06dbeSstefano_zampini BDdelta_DN ctx; 7333ca39a21SBarry Smith MatSolverType solver; 73425d06dbeSstefano_zampini const char *prefix; 73525d06dbeSstefano_zampini 7369566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat_ctx->B_Ddelta, &ctx)); 7379566063dSJacob Faibussowitsch PetscCall(KSPSetType(ctx->kBD, KSPPREONLY)); 7389566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ctx->kBD, &mpc)); 7399566063dSJacob Faibussowitsch PetscCall(KSPGetPC(pcbddc->ksp_D, &pc)); 7409566063dSJacob Faibussowitsch PetscCall(PCSetType(mpc, PCLU)); 741835f2295SStefano Zampini PetscCall(PCFactorGetMatSolverType(pc, &solver)); 7421baa6e33SBarry Smith if (solver) PetscCall(PCFactorSetMatSolverType(mpc, solver)); 7439566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(fetimat, &prefix)); 7449566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ctx->kBD, prefix)); 7459566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ctx->kBD, "bddelta_")); 7469566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ctx->kBD)); 74725d06dbeSstefano_zampini } 74825d06dbeSstefano_zampini 749e1214c54Sstefano_zampini if (mat_ctx->l2g_lambda_only) { 7509566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only)); 751e1214c54Sstefano_zampini fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only; 752e1214c54Sstefano_zampini } else { 7539566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda)); 754674ae819SStefano Zampini fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda; 755e1214c54Sstefano_zampini } 756f28b6018SStefano Zampini /* Dirichlet preconditioner */ 7579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_lumped", &lumped, NULL)); 758f28b6018SStefano Zampini if (!lumped) { 7599feb908dSstefano_zampini IS iV; 7609c2d02cdSstefano_zampini PetscBool discrete_harmonic = PETSC_FALSE; 7619c2d02cdSstefano_zampini 7629566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)fetidppc_ctx->pc, "__KSPFETIDP_iV", (PetscObject *)&iV)); 76348a46eb9SPierre Jolivet if (iV) PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_discrete_harmonic", &discrete_harmonic, NULL)); 7649c2d02cdSstefano_zampini if (discrete_harmonic) { 7659c2d02cdSstefano_zampini KSP sksp; 7669c2d02cdSstefano_zampini PC pc; 767111315fdSstefano_zampini PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs; 7689c2d02cdSstefano_zampini Mat A_II, A_IB, A_BI; 769111315fdSstefano_zampini IS iP = NULL; 770111315fdSstefano_zampini PetscBool isshell, reuse = PETSC_FALSE; 7719c2d02cdSstefano_zampini KSPType ksptype; 7729c2d02cdSstefano_zampini const char *prefix; 7739c2d02cdSstefano_zampini 7749c2d02cdSstefano_zampini /* 7759c2d02cdSstefano_zampini We constructs a Schur complement for 7769c2d02cdSstefano_zampini 7779c2d02cdSstefano_zampini | A_II A_ID | 7789c2d02cdSstefano_zampini | A_DI A_DD | 7799c2d02cdSstefano_zampini 7809c2d02cdSstefano_zampini instead of 7819c2d02cdSstefano_zampini 7829c2d02cdSstefano_zampini | A_II B^t_II A_ID | 7839c2d02cdSstefano_zampini | B_II -C_II B_ID | 7849c2d02cdSstefano_zampini | A_DI B^t_ID A_DD | 7859c2d02cdSstefano_zampini 7869c2d02cdSstefano_zampini */ 787111315fdSstefano_zampini if (sub_schurs && sub_schurs->reuse_solver) { 7889566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)sub_schurs->A, "__KSPFETIDP_iP", (PetscObject *)&iP)); 789111315fdSstefano_zampini if (iP) reuse = PETSC_TRUE; 790111315fdSstefano_zampini } 791111315fdSstefano_zampini if (!reuse) { 792111315fdSstefano_zampini IS aB; 793111315fdSstefano_zampini PetscInt nb; 7949566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(pcis->is_B_local, &nb)); 7959566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II), nb, 0, 1, &aB)); 7969566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pcis->A_II, iV, iV, MAT_INITIAL_MATRIX, &A_II)); 7979566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pcis->A_IB, iV, aB, MAT_INITIAL_MATRIX, &A_IB)); 7989566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pcis->A_BI, aB, iV, MAT_INITIAL_MATRIX, &A_BI)); 7999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&aB)); 800111315fdSstefano_zampini } else { 8019566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_I_local, pcis->is_B_local, MAT_INITIAL_MATRIX, &A_IB)); 8029566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_B_local, pcis->is_I_local, MAT_INITIAL_MATRIX, &A_BI)); 8039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pcis->A_II)); 804111315fdSstefano_zampini A_II = pcis->A_II; 805111315fdSstefano_zampini } 8069566063dSJacob Faibussowitsch PetscCall(MatCreateSchurComplement(A_II, A_II, A_IB, A_BI, pcis->A_BB, &fetidppc_ctx->S_j)); 8079c2d02cdSstefano_zampini 8089c2d02cdSstefano_zampini /* propagate settings of solver */ 8099566063dSJacob Faibussowitsch PetscCall(MatSchurComplementGetKSP(fetidppc_ctx->S_j, &sksp)); 8109566063dSJacob Faibussowitsch PetscCall(KSPGetType(pcis->ksp_D, &ksptype)); 8119566063dSJacob Faibussowitsch PetscCall(KSPSetType(sksp, ksptype)); 8129566063dSJacob Faibussowitsch PetscCall(KSPGetPC(pcis->ksp_D, &pc)); 8139566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &isshell)); 8145334bea6Sstefano_zampini if (!isshell) { 8153ca39a21SBarry Smith MatSolverType solver; 8165334bea6Sstefano_zampini PCType pctype; 8175334bea6Sstefano_zampini 8189566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &pctype)); 819835f2295SStefano Zampini PetscCall(PCFactorGetMatSolverType(pc, &solver)); 8209566063dSJacob Faibussowitsch PetscCall(KSPGetPC(sksp, &pc)); 8219566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, pctype)); 8221baa6e33SBarry Smith if (solver) PetscCall(PCFactorSetMatSolverType(pc, solver)); 8235334bea6Sstefano_zampini } else { 8249566063dSJacob Faibussowitsch PetscCall(KSPGetPC(sksp, &pc)); 8259566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCLU)); 8265334bea6Sstefano_zampini } 8279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_II)); 8289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_IB)); 8299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_BI)); 8309566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(fetimat, &prefix)); 8319566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(sksp, prefix)); 8329566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(sksp, "harmonic_")); 8339566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(sksp)); 834111315fdSstefano_zampini if (reuse) { 8359566063dSJacob Faibussowitsch PetscCall(KSPSetPC(sksp, sub_schurs->reuse_solver->interior_solver)); 8369566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver, (PetscObject)sksp, 0)); 837111315fdSstefano_zampini } 8389c2d02cdSstefano_zampini } else { /* default Dirichlet preconditioner is pde-harmonic */ 8399566063dSJacob Faibussowitsch PetscCall(MatCreateSchurComplement(pcis->A_II, pcis->A_II, pcis->A_IB, pcis->A_BI, pcis->A_BB, &fetidppc_ctx->S_j)); 8409566063dSJacob Faibussowitsch PetscCall(MatSchurComplementSetKSP(fetidppc_ctx->S_j, pcis->ksp_D)); 8419c2d02cdSstefano_zampini } 842f28b6018SStefano Zampini } else { 8439566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pcis->A_BB)); 844f28b6018SStefano Zampini fetidppc_ctx->S_j = pcis->A_BB; 845f28b6018SStefano Zampini } 846af140850Sstefano_zampini /* saddle-point */ 847af140850Sstefano_zampini if (mat_ctx->xPg) { 8489566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->xPg)); 849af140850Sstefano_zampini fetidppc_ctx->xPg = mat_ctx->xPg; 8509566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat_ctx->yPg)); 851af140850Sstefano_zampini fetidppc_ctx->yPg = mat_ctx->yPg; 852af140850Sstefano_zampini } 8533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 854674ae819SStefano Zampini } 855674ae819SStefano Zampini 85666976f2fSJacob Faibussowitsch static PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans) 857d71ae5a4SJacob Faibussowitsch { 858674ae819SStefano Zampini FETIDPMat_ctx mat_ctx; 859617d11aeSStefano Zampini PC_BDDC *pcbddc; 860674ae819SStefano Zampini PC_IS *pcis; 861674ae819SStefano Zampini 862674ae819SStefano Zampini PetscFunctionBegin; 8639566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(fetimat, &mat_ctx)); 864674ae819SStefano Zampini pcis = (PC_IS *)mat_ctx->pc->data; 865617d11aeSStefano Zampini pcbddc = (PC_BDDC *)mat_ctx->pc->data; 866674ae819SStefano Zampini /* Application of B_delta^T */ 8679566063dSJacob Faibussowitsch PetscCall(VecSet(pcis->vec1_B, 0.)); 8689566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE)); 8699566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE)); 8709566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat_ctx->B_delta, mat_ctx->lambda_local, pcis->vec1_B)); 871af140850Sstefano_zampini 872af140850Sstefano_zampini /* Add contribution from saddle point */ 873af140850Sstefano_zampini if (mat_ctx->l2g_p) { 8749566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE)); 8759566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE)); 876af140850Sstefano_zampini if (pcbddc->switch_static) { 877ef425aeaSstefano_zampini if (trans) { 8789566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat_ctx->B_BI, mat_ctx->vP, pcis->vec1_D)); 879ef425aeaSstefano_zampini } else { 8809566063dSJacob Faibussowitsch PetscCall(MatMult(mat_ctx->Bt_BI, mat_ctx->vP, pcis->vec1_D)); 881af140850Sstefano_zampini } 882ef425aeaSstefano_zampini } 883ef425aeaSstefano_zampini if (trans) { 8849566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(mat_ctx->B_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B)); 885ef425aeaSstefano_zampini } else { 8869566063dSJacob Faibussowitsch PetscCall(MatMultAdd(mat_ctx->Bt_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B)); 887ef425aeaSstefano_zampini } 888af140850Sstefano_zampini } else { 8891baa6e33SBarry Smith if (pcbddc->switch_static) PetscCall(VecSet(pcis->vec1_D, 0.0)); 890af140850Sstefano_zampini } 891af140850Sstefano_zampini /* Application of \widetilde{S}^-1 */ 8929566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n)); 8939566063dSJacob Faibussowitsch PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc, trans)); 8949566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n)); 8959566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 896674ae819SStefano Zampini /* Application of B_delta */ 8979566063dSJacob Faibussowitsch PetscCall(MatMult(mat_ctx->B_delta, pcis->vec1_B, mat_ctx->lambda_local)); 898af140850Sstefano_zampini /* Contribution from boundary pressures */ 899af140850Sstefano_zampini if (mat_ctx->C) { 900af140850Sstefano_zampini const PetscScalar *lx; 901af140850Sstefano_zampini PetscScalar *ly; 902af140850Sstefano_zampini 90315579a77SStefano Zampini /* pressure ordered first in the local part of x and y */ 9049566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &lx)); 9059566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &ly)); 9069566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(mat_ctx->xPg, lx)); 9079566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(mat_ctx->yPg, ly)); 908ef425aeaSstefano_zampini if (trans) { 9099566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg)); 910ef425aeaSstefano_zampini } else { 9119566063dSJacob Faibussowitsch PetscCall(MatMult(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg)); 912ef425aeaSstefano_zampini } 9139566063dSJacob Faibussowitsch PetscCall(VecResetArray(mat_ctx->xPg)); 9149566063dSJacob Faibussowitsch PetscCall(VecResetArray(mat_ctx->yPg)); 9159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &lx)); 9169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &ly)); 917af140850Sstefano_zampini } 918af140850Sstefano_zampini /* Add contribution from saddle point */ 919af140850Sstefano_zampini if (mat_ctx->l2g_p) { 92026699680SStefano Zampini PetscCall(VecISSet(pcis->vec1_B, mat_ctx->lP_B, 0)); 921ef425aeaSstefano_zampini if (trans) { 9229566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat_ctx->Bt_BB, pcis->vec1_B, mat_ctx->vP)); 923ef425aeaSstefano_zampini } else { 9249566063dSJacob Faibussowitsch PetscCall(MatMult(mat_ctx->B_BB, pcis->vec1_B, mat_ctx->vP)); 925ef425aeaSstefano_zampini } 926af140850Sstefano_zampini if (pcbddc->switch_static) { 92726699680SStefano Zampini PetscCall(VecISSet(pcis->vec1_D, mat_ctx->lP_I, 0)); 928ef425aeaSstefano_zampini if (trans) { 9299566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(mat_ctx->Bt_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP)); 930ef425aeaSstefano_zampini } else { 9319566063dSJacob Faibussowitsch PetscCall(MatMultAdd(mat_ctx->B_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP)); 932af140850Sstefano_zampini } 933ef425aeaSstefano_zampini } 9349566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD)); 9359566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD)); 936af140850Sstefano_zampini } 9379566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD)); 9389566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD)); 9393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 940674ae819SStefano Zampini } 941674ae819SStefano Zampini 942d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) 943d71ae5a4SJacob Faibussowitsch { 944edf7251bSStefano Zampini PetscFunctionBegin; 9459566063dSJacob Faibussowitsch PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_FALSE)); 9463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 947ef425aeaSstefano_zampini } 948ef425aeaSstefano_zampini 949d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) 950d71ae5a4SJacob Faibussowitsch { 951ef425aeaSstefano_zampini PetscFunctionBegin; 9529566063dSJacob Faibussowitsch PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_TRUE)); 9533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 954edf7251bSStefano Zampini } 955edf7251bSStefano Zampini 95666976f2fSJacob Faibussowitsch static PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans) 957d71ae5a4SJacob Faibussowitsch { 958674ae819SStefano Zampini FETIDPPC_ctx pc_ctx; 959674ae819SStefano Zampini PC_IS *pcis; 960674ae819SStefano Zampini 961674ae819SStefano Zampini PetscFunctionBegin; 9629566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(fetipc, &pc_ctx)); 963674ae819SStefano Zampini pcis = (PC_IS *)pc_ctx->pc->data; 964674ae819SStefano Zampini /* Application of B_Ddelta^T */ 9659566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE)); 9669566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE)); 9679566063dSJacob Faibussowitsch PetscCall(VecSet(pcis->vec2_B, 0.0)); 9689566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc_ctx->B_Ddelta, pc_ctx->lambda_local, pcis->vec2_B)); 969ed6c3d69SStefano Zampini /* Application of local Schur complement */ 97016597391Sstefano_zampini if (trans) { 9719566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B)); 97216597391Sstefano_zampini } else { 9739566063dSJacob Faibussowitsch PetscCall(MatMult(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B)); 97416597391Sstefano_zampini } 975edf7251bSStefano Zampini /* Application of B_Ddelta */ 9769566063dSJacob Faibussowitsch PetscCall(MatMult(pc_ctx->B_Ddelta, pcis->vec1_B, pc_ctx->lambda_local)); 9779566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 9789566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD)); 9799566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD)); 9803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 981edf7251bSStefano Zampini } 982edf7251bSStefano Zampini 983d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y) 984d71ae5a4SJacob Faibussowitsch { 985edf7251bSStefano Zampini PetscFunctionBegin; 9869566063dSJacob Faibussowitsch PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_FALSE)); 9873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 98816597391Sstefano_zampini } 98916597391Sstefano_zampini 990d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y) 991d71ae5a4SJacob Faibussowitsch { 99216597391Sstefano_zampini PetscFunctionBegin; 9939566063dSJacob Faibussowitsch PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_TRUE)); 9943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 995674ae819SStefano Zampini } 996c45b8d2dSstefano_zampini 997d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer) 998d71ae5a4SJacob Faibussowitsch { 999c45b8d2dSstefano_zampini FETIDPPC_ctx pc_ctx; 10009f196a02SMartin Diehl PetscBool isascii; 1001c45b8d2dSstefano_zampini PetscViewer sviewer; 1002c45b8d2dSstefano_zampini 1003c45b8d2dSstefano_zampini PetscFunctionBegin; 10049f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 10059f196a02SMartin Diehl if (isascii) { 1006c45b8d2dSstefano_zampini PetscMPIInt rank; 100725d06dbeSstefano_zampini PetscBool isschur, isshell; 1008c45b8d2dSstefano_zampini 10099566063dSJacob Faibussowitsch PetscCall(PCShellGetContext(pc, &pc_ctx)); 10109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 10119566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->S_j, MATSCHURCOMPLEMENT, &isschur)); 1012c45b8d2dSstefano_zampini if (isschur) { 10139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Dirichlet preconditioner (just from rank 0)\n")); 1014c45b8d2dSstefano_zampini } else { 10159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Lumped preconditioner (just from rank 0)\n")); 1016c45b8d2dSstefano_zampini } 10179566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer)); 1018dd400576SPatrick Sanan if (rank == 0) { 10199566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO)); 10209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(sviewer)); 10219566063dSJacob Faibussowitsch PetscCall(MatView(pc_ctx->S_j, sviewer)); 10229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(sviewer)); 10239566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(sviewer)); 1024c45b8d2dSstefano_zampini } 10259566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer)); 10269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta, MATSHELL, &isshell)); 102725d06dbeSstefano_zampini if (isshell) { 102825d06dbeSstefano_zampini BDdelta_DN ctx; 10299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n")); 10309566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(pc_ctx->B_Ddelta, &ctx)); 10319566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer)); 1032dd400576SPatrick Sanan if (rank == 0) { 103315579a77SStefano Zampini PetscInt tl; 103415579a77SStefano Zampini 10359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetTab(sviewer, &tl)); 10369566063dSJacob Faibussowitsch PetscCall(PetscObjectSetTabLevel((PetscObject)ctx->kBD, tl)); 10379566063dSJacob Faibussowitsch PetscCall(KSPView(ctx->kBD, sviewer)); 10389566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO)); 10399566063dSJacob Faibussowitsch PetscCall(MatView(ctx->BD, sviewer)); 10409566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(sviewer)); 104125d06dbeSstefano_zampini } 10429566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer)); 1043e5afcf28SBarry Smith } 10449566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1045c45b8d2dSstefano_zampini } 10463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1047c45b8d2dSstefano_zampini } 1048