132fe681dSStefano Zampini #include <petsc/private/pcmgimpl.h> 232fe681dSStefano Zampini #include <petsc/private/pcbddcimpl.h> 332fe681dSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h> 432fe681dSStefano Zampini 5d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMGGDSWSetUp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat A, PetscInt *ns, Mat **sA_IG_n, KSP **sksp_n, IS **sI_n, IS **sG_n, Mat **sGf_n, IS **sGi_n, IS **sGiM_n) 6d71ae5a4SJacob Faibussowitsch { 732fe681dSStefano Zampini KSP *sksp; 832fe681dSStefano Zampini PC pcbddc = NULL, smoothpc; 932fe681dSStefano Zampini PC_BDDC *ipcbddc; 1032fe681dSStefano Zampini PC_IS *ipcis; 1132fe681dSStefano Zampini Mat *sA_IG, *sGf, cmat, lA; 1232fe681dSStefano Zampini ISLocalToGlobalMapping l2g; 1332fe681dSStefano Zampini IS *sI, *sG, *sGi, *sGiM, cref; 1432fe681dSStefano Zampini PCBDDCSubSchurs sub_schurs = NULL; 1532fe681dSStefano Zampini PCBDDCGraph graph; 1632fe681dSStefano Zampini const char *prefix; 1732fe681dSStefano Zampini const PetscScalar *tdata; 1832fe681dSStefano Zampini PetscScalar *data, *cdata; 1932fe681dSStefano Zampini PetscReal tol = 0.0, otol; 2032fe681dSStefano Zampini const PetscInt *ia, *ja; 2132fe681dSStefano Zampini PetscInt *ccii, *cridx; 221690c2aeSBarry Smith PetscInt i, j, ngct, ng, dbg = 0, odbg, minmax[2] = {0, PETSC_INT_MAX}, ominmax[2], vsize; 2332fe681dSStefano Zampini PetscBool flg, userdefined = PETSC_TRUE, reuse_solver = PETSC_TRUE, reduced = PETSC_FALSE; 2432fe681dSStefano Zampini 2532fe681dSStefano Zampini PetscFunctionBegin; 2632fe681dSStefano Zampini PetscCall(MatGetBlockSize(A, &vsize)); 2732fe681dSStefano Zampini PetscCall(KSPGetOptionsPrefix(smooth, &prefix)); 2832fe681dSStefano Zampini PetscOptionsBegin(PetscObjectComm((PetscObject)smooth), prefix, "GDSW options", "PC"); 2932fe681dSStefano Zampini PetscCall(PetscOptionsReal("-gdsw_tolerance", "Tolerance for eigenvalue problem", NULL, tol, &tol, NULL)); 3032fe681dSStefano Zampini PetscCall(PetscOptionsBool("-gdsw_userdefined", "Use user-defined functions in addition to those adaptively generated", NULL, userdefined, &userdefined, NULL)); 3132fe681dSStefano Zampini PetscCall(PetscOptionsIntArray("-gdsw_minmax", "Minimum and maximum number of basis functions per connected component for adaptive GDSW", NULL, minmax, (i = 2, &i), NULL)); 3232fe681dSStefano Zampini PetscCall(PetscOptionsInt("-gdsw_vertex_size", "Connected components smaller or equal to vertex size will be considered as vertices", NULL, vsize, &vsize, NULL)); 3332fe681dSStefano Zampini PetscCall(PetscOptionsBool("-gdsw_reuse", "Reuse interior solver from Schur complement computations", NULL, reuse_solver, &reuse_solver, NULL)); 3432fe681dSStefano Zampini PetscCall(PetscOptionsBool("-gdsw_reduced", "Reduced GDSW", NULL, reduced, &reduced, NULL)); 3532fe681dSStefano Zampini PetscCall(PetscOptionsInt("-gdsw_debug", "Debug output", NULL, dbg, &dbg, NULL)); 3632fe681dSStefano Zampini PetscOptionsEnd(); 3732fe681dSStefano Zampini 3832fe681dSStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)A, MATIS, &flg)); 3932fe681dSStefano Zampini if (!flg) { 4032fe681dSStefano Zampini MatNullSpace nnsp; 4132fe681dSStefano Zampini 4232fe681dSStefano Zampini PetscCall(MatGetNearNullSpace(A, &nnsp)); 433ba16761SJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)nnsp)); 4432fe681dSStefano Zampini PetscCall(MatConvert(A, MATIS, MAT_INITIAL_MATRIX, &A)); 4532fe681dSStefano Zampini PetscCall(MatSetNearNullSpace(A, nnsp)); 4632fe681dSStefano Zampini PetscCall(MatNullSpaceDestroy(&nnsp)); 4732fe681dSStefano Zampini } else PetscCall(PetscObjectReference((PetscObject)A)); 4832fe681dSStefano Zampini 4932fe681dSStefano Zampini /* TODO Multi sub */ 5032fe681dSStefano Zampini *ns = 1; 5132fe681dSStefano Zampini PetscCall(PetscMalloc1(*ns, &sA_IG)); 5232fe681dSStefano Zampini PetscCall(PetscMalloc1(*ns, &sksp)); 5332fe681dSStefano Zampini PetscCall(PetscMalloc1(*ns, &sI)); 5432fe681dSStefano Zampini PetscCall(PetscMalloc1(*ns, &sG)); 5532fe681dSStefano Zampini PetscCall(PetscMalloc1(*ns, &sGf)); 5632fe681dSStefano Zampini PetscCall(PetscMalloc1(*ns, &sGi)); 5732fe681dSStefano Zampini PetscCall(PetscMalloc1(*ns, &sGiM)); 5832fe681dSStefano Zampini *sA_IG_n = sA_IG; 5932fe681dSStefano Zampini *sksp_n = sksp; 6032fe681dSStefano Zampini *sI_n = sI; 6132fe681dSStefano Zampini *sG_n = sG; 6232fe681dSStefano Zampini *sGf_n = sGf; 6332fe681dSStefano Zampini *sGi_n = sGi; 6432fe681dSStefano Zampini *sGiM_n = sGiM; 6532fe681dSStefano Zampini 6632fe681dSStefano Zampini /* submatrices and solvers */ 6732fe681dSStefano Zampini PetscCall(KSPGetPC(smooth, &smoothpc)); 6832fe681dSStefano Zampini PetscCall(PetscObjectTypeCompareAny((PetscObject)smoothpc, &flg, PCBDDC, "")); 6932fe681dSStefano Zampini if (!flg) { 7032fe681dSStefano Zampini Mat smoothA; 7132fe681dSStefano Zampini 7232fe681dSStefano Zampini PetscCall(PCGetOperators(smoothpc, &smoothA, NULL)); 7332fe681dSStefano Zampini PetscCall(PCCreate(PetscObjectComm((PetscObject)A), &pcbddc)); 7432fe681dSStefano Zampini PetscCall(PCSetType(pcbddc, PCBDDC)); 7532fe681dSStefano Zampini PetscCall(PCSetOperators(pcbddc, smoothA, A)); 7632fe681dSStefano Zampini PetscCall(PCISSetUp(pcbddc, PETSC_TRUE, PETSC_FALSE)); 7732fe681dSStefano Zampini } else { 7832fe681dSStefano Zampini PetscCall(PetscObjectReference((PetscObject)smoothpc)); 7932fe681dSStefano Zampini pcbddc = smoothpc; 8032fe681dSStefano Zampini } 8132fe681dSStefano Zampini ipcis = (PC_IS *)pcbddc->data; 8232fe681dSStefano Zampini ipcbddc = (PC_BDDC *)pcbddc->data; 8332fe681dSStefano Zampini PetscCall(PetscObjectReference((PetscObject)ipcis->A_IB)); 8432fe681dSStefano Zampini PetscCall(PetscObjectReference((PetscObject)ipcis->is_I_global)); 8532fe681dSStefano Zampini PetscCall(PetscObjectReference((PetscObject)ipcis->is_B_global)); 8632fe681dSStefano Zampini sA_IG[0] = ipcis->A_IB; 8732fe681dSStefano Zampini sI[0] = ipcis->is_I_global; 8832fe681dSStefano Zampini sG[0] = ipcis->is_B_global; 8932fe681dSStefano Zampini 9032fe681dSStefano Zampini PetscCall(KSPCreate(PetscObjectComm((PetscObject)ipcis->A_II), &sksp[0])); 913821be0aSBarry Smith PetscCall(KSPSetNestLevel(sksp[0], pc->kspnestlevel)); 9232fe681dSStefano Zampini PetscCall(KSPSetOperators(sksp[0], ipcis->A_II, ipcis->pA_II)); 9332fe681dSStefano Zampini PetscCall(KSPSetOptionsPrefix(sksp[0], prefix)); 9432fe681dSStefano Zampini PetscCall(KSPAppendOptionsPrefix(sksp[0], "gdsw_")); 9532fe681dSStefano Zampini PetscCall(KSPSetFromOptions(sksp[0])); 9632fe681dSStefano Zampini 9732fe681dSStefano Zampini /* analyze interface */ 9832fe681dSStefano Zampini PetscCall(MatISGetLocalMat(A, &lA)); 9932fe681dSStefano Zampini graph = ipcbddc->mat_graph; 10032fe681dSStefano Zampini if (!flg) { 10132fe681dSStefano Zampini PetscInt N; 10232fe681dSStefano Zampini 10332fe681dSStefano Zampini PetscCall(MatISGetLocalToGlobalMapping(A, &l2g, NULL)); 10432fe681dSStefano Zampini PetscCall(MatGetSize(A, &N, NULL)); 1051690c2aeSBarry Smith PetscCall(PCBDDCGraphInit(graph, l2g, N, PETSC_INT_MAX)); 10632fe681dSStefano Zampini PetscCall(MatGetRowIJ(lA, 0, PETSC_TRUE, PETSC_FALSE, &graph->nvtxs_csr, (const PetscInt **)&graph->xadj, (const PetscInt **)&graph->adjncy, &flg)); 10732fe681dSStefano Zampini PetscCall(PCBDDCGraphSetUp(graph, vsize, NULL, NULL, 0, NULL, NULL)); 10832fe681dSStefano Zampini PetscCall(MatRestoreRowIJ(lA, 0, PETSC_TRUE, PETSC_FALSE, &graph->nvtxs_csr, (const PetscInt **)&graph->xadj, (const PetscInt **)&graph->adjncy, &flg)); 10932fe681dSStefano Zampini PetscCall(PCBDDCGraphComputeConnectedComponents(graph)); 11032fe681dSStefano Zampini } 11132fe681dSStefano Zampini l2g = graph->l2gmap; 11232fe681dSStefano Zampini if (reduced) { 11332fe681dSStefano Zampini PetscContainer gcand; 11432fe681dSStefano Zampini PCBDDCGraphCandidates cand; 11532fe681dSStefano Zampini PetscErrorCode (*rgdsw)(DM, PetscInt *, IS **); 11632fe681dSStefano Zampini 11732fe681dSStefano Zampini PetscCall(PetscObjectQueryFunction((PetscObject)dm, "DMComputeLocalRGDSWSets", &rgdsw)); 11832fe681dSStefano Zampini PetscCheck(rgdsw, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not supported"); 11932fe681dSStefano Zampini PetscCall(PetscNew(&cand)); 12032fe681dSStefano Zampini PetscCall((*rgdsw)(dm, &cand->nfc, &cand->Faces)); 12132fe681dSStefano Zampini /* filter interior (if any) and guarantee IS are ordered by global numbering */ 12232fe681dSStefano Zampini for (i = 0; i < cand->nfc; i++) { 12332fe681dSStefano Zampini IS is, is2; 12432fe681dSStefano Zampini 12532fe681dSStefano Zampini PetscCall(ISLocalToGlobalMappingApplyIS(l2g, cand->Faces[i], &is)); 12632fe681dSStefano Zampini PetscCall(ISDestroy(&cand->Faces[i])); 12732fe681dSStefano Zampini PetscCall(ISSort(is)); 12832fe681dSStefano Zampini PetscCall(ISGlobalToLocalMappingApplyIS(l2g, IS_GTOLM_DROP, is, &is2)); 12932fe681dSStefano Zampini PetscCall(ISDestroy(&is)); 13032fe681dSStefano Zampini PetscCall(ISGlobalToLocalMappingApplyIS(ipcis->BtoNmap, IS_GTOLM_DROP, is2, &is)); 13132fe681dSStefano Zampini PetscCall(ISDestroy(&is2)); 13232fe681dSStefano Zampini PetscCall(ISLocalToGlobalMappingApplyIS(ipcis->BtoNmap, is, &cand->Faces[i])); 13332fe681dSStefano Zampini PetscCall(ISDestroy(&is)); 13432fe681dSStefano Zampini } 13532fe681dSStefano Zampini PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &gcand)); 13632fe681dSStefano Zampini PetscCall(PetscContainerSetPointer(gcand, cand)); 137*49abdd8aSBarry Smith PetscCall(PetscContainerSetCtxDestroy(gcand, PCBDDCDestroyGraphCandidatesIS)); 13832fe681dSStefano Zampini PetscCall(PetscObjectCompose((PetscObject)l2g, "_PCBDDCGraphCandidatesIS", (PetscObject)gcand)); 13932fe681dSStefano Zampini PetscCall(PetscContainerDestroy(&gcand)); 14032fe681dSStefano Zampini } 14132fe681dSStefano Zampini 14232fe681dSStefano Zampini /* interface functions */ 14332fe681dSStefano Zampini otol = ipcbddc->adaptive_threshold[1]; 14432fe681dSStefano Zampini odbg = ipcbddc->dbg_flag; 14532fe681dSStefano Zampini ominmax[0] = ipcbddc->adaptive_nmin; 14632fe681dSStefano Zampini ominmax[1] = ipcbddc->adaptive_nmax; 14732fe681dSStefano Zampini ipcbddc->adaptive_threshold[1] = tol; 14832fe681dSStefano Zampini ipcbddc->dbg_flag = dbg; 14932fe681dSStefano Zampini ipcbddc->adaptive_nmin = minmax[0]; 15032fe681dSStefano Zampini ipcbddc->adaptive_nmax = minmax[1]; 15132fe681dSStefano Zampini if (tol != 0.0) { /* adaptive */ 15232fe681dSStefano Zampini Mat lS; 15332fe681dSStefano Zampini 15432fe681dSStefano Zampini PetscCall(MatCreateSchurComplement(ipcis->A_II, ipcis->pA_II, ipcis->A_IB, ipcis->A_BI, ipcis->A_BB, &lS)); 15532fe681dSStefano Zampini PetscCall(KSPGetOptionsPrefix(sksp[0], &prefix)); 15632fe681dSStefano Zampini PetscCall(PCBDDCSubSchursCreate(&sub_schurs)); 15732fe681dSStefano Zampini PetscCall(PCBDDCSubSchursInit(sub_schurs, prefix, ipcis->is_I_local, ipcis->is_B_local, graph, ipcis->BtoNmap, PETSC_FALSE, PETSC_TRUE)); 15832fe681dSStefano Zampini if (userdefined) PetscCall(PCBDDCComputeFakeChange(pcbddc, PETSC_FALSE, graph, NULL, &cmat, &cref, NULL, &flg)); 1599371c9d4SSatish Balay else { 1609371c9d4SSatish Balay cmat = NULL; 1619371c9d4SSatish Balay cref = NULL; 1629371c9d4SSatish Balay } 16332fe681dSStefano Zampini PetscCall(PCBDDCSubSchursSetUp(sub_schurs, lA, lS, PETSC_TRUE, NULL, NULL, -1, NULL, PETSC_TRUE, reuse_solver, PETSC_FALSE, 0, NULL, NULL, cmat, cref)); 16432fe681dSStefano Zampini PetscCall(MatDestroy(&lS)); 16532fe681dSStefano Zampini PetscCall(MatDestroy(&cmat)); 16632fe681dSStefano Zampini PetscCall(ISDestroy(&cref)); 16732fe681dSStefano Zampini if (sub_schurs->reuse_solver) { 16832fe681dSStefano Zampini PetscCall(KSPSetPC(sksp[0], sub_schurs->reuse_solver->interior_solver)); 16932fe681dSStefano Zampini PetscCall(PCDestroy(&sub_schurs->reuse_solver->interior_solver)); 17032fe681dSStefano Zampini sub_schurs->reuse_solver = NULL; 17132fe681dSStefano Zampini } 17232fe681dSStefano Zampini } 17332fe681dSStefano Zampini PetscCall(PCBDDCComputeFakeChange(pcbddc, PETSC_TRUE, graph, sub_schurs, &cmat, &cref, &sGiM[0], NULL)); 17432fe681dSStefano Zampini PetscCall(PCBDDCSubSchursDestroy(&sub_schurs)); 17532fe681dSStefano Zampini ipcbddc->adaptive_threshold[1] = otol; 17632fe681dSStefano Zampini ipcbddc->dbg_flag = odbg; 17732fe681dSStefano Zampini ipcbddc->adaptive_nmin = ominmax[0]; 17832fe681dSStefano Zampini ipcbddc->adaptive_nmax = ominmax[1]; 17932fe681dSStefano Zampini 18032fe681dSStefano Zampini PetscCall(ISLocalToGlobalMappingApplyIS(l2g, cref, &sGi[0])); 18132fe681dSStefano Zampini PetscCall(ISDestroy(&cref)); 18232fe681dSStefano Zampini 18332fe681dSStefano Zampini PetscCall(MatSeqAIJGetArrayRead(cmat, &tdata)); 18432fe681dSStefano Zampini PetscCall(MatGetRowIJ(cmat, 0, PETSC_FALSE, PETSC_FALSE, &ngct, &ia, &ja, &flg)); 18532fe681dSStefano Zampini PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in MatGetRowIJ"); 18632fe681dSStefano Zampini 18732fe681dSStefano Zampini PetscCall(PetscMalloc1(ngct + 1, &ccii)); 18832fe681dSStefano Zampini PetscCall(PetscMalloc1(ia[ngct], &cridx)); 18932fe681dSStefano Zampini PetscCall(PetscMalloc1(ia[ngct], &cdata)); 19032fe681dSStefano Zampini 19132fe681dSStefano Zampini PetscCall(PetscArraycpy(ccii, ia, ngct + 1)); 19232fe681dSStefano Zampini PetscCall(PetscArraycpy(cdata, tdata, ia[ngct])); 19332fe681dSStefano Zampini PetscCall(ISGlobalToLocalMappingApply(ipcis->BtoNmap, IS_GTOLM_DROP, ia[ngct], ja, &i, cridx)); 19432fe681dSStefano Zampini PetscCheck(i == ia[ngct], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in G2L"); 19532fe681dSStefano Zampini 19632fe681dSStefano Zampini PetscCall(MatRestoreRowIJ(cmat, 0, PETSC_FALSE, PETSC_FALSE, &i, &ia, &ja, &flg)); 19732fe681dSStefano Zampini PetscCall(MatSeqAIJRestoreArrayRead(cmat, &tdata)); 19832fe681dSStefano Zampini PetscCall(MatDestroy(&cmat)); 19932fe681dSStefano Zampini 20032fe681dSStefano Zampini /* populate dense matrix */ 20132fe681dSStefano Zampini PetscCall(ISGetLocalSize(sG[0], &ng)); 20232fe681dSStefano Zampini PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, ng, ngct, NULL, &sGf[0])); 20332fe681dSStefano Zampini PetscCall(MatDenseGetArrayWrite(sGf[0], &data)); 20432fe681dSStefano Zampini for (i = 0; i < ngct; i++) 2059371c9d4SSatish Balay for (j = ccii[i]; j < ccii[i + 1]; j++) data[ng * i + cridx[j]] = cdata[j]; 20632fe681dSStefano Zampini PetscCall(MatDenseRestoreArrayWrite(sGf[0], &data)); 20732fe681dSStefano Zampini 20832fe681dSStefano Zampini PetscCall(PetscFree(cdata)); 20932fe681dSStefano Zampini PetscCall(PetscFree(ccii)); 21032fe681dSStefano Zampini PetscCall(PetscFree(cridx)); 21132fe681dSStefano Zampini PetscCall(PCDestroy(&pcbddc)); 21232fe681dSStefano Zampini PetscCall(MatDestroy(&A)); 2133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21432fe681dSStefano Zampini } 21532fe681dSStefano Zampini 216d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGDSWCreateCoarseSpace_Private(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat guess, Mat *cspace) 217d71ae5a4SJacob Faibussowitsch { 21832fe681dSStefano Zampini KSP *sksp; 21932fe681dSStefano Zampini Mat A, *sA_IG, *sGf, preallocator; 22032fe681dSStefano Zampini IS Gidx, GidxMult, cG; 22132fe681dSStefano Zampini IS *sI, *sG, *sGi, *sGiM; 22232fe681dSStefano Zampini const PetscInt *cidx; 22332fe681dSStefano Zampini PetscInt NG, ns, n, i, c, rbs, cbs[2]; 22432fe681dSStefano Zampini PetscBool flg; 22532fe681dSStefano Zampini MatType ptype; 22632fe681dSStefano Zampini 22732fe681dSStefano Zampini PetscFunctionBegin; 22832fe681dSStefano Zampini *cspace = NULL; 2293ba16761SJacob Faibussowitsch if (!l) PetscFunctionReturn(PETSC_SUCCESS); 23032fe681dSStefano Zampini if (pc->useAmat) { 23132fe681dSStefano Zampini PetscCall(KSPGetOperatorsSet(smooth, &flg, NULL)); 23232fe681dSStefano Zampini PetscCheck(flg, PetscObjectComm((PetscObject)smooth), PETSC_ERR_ORDER, "Amat not set"); 23332fe681dSStefano Zampini PetscCall(KSPGetOperators(smooth, &A, NULL)); 23432fe681dSStefano Zampini } else { 23532fe681dSStefano Zampini PetscCall(KSPGetOperatorsSet(smooth, NULL, &flg)); 23632fe681dSStefano Zampini PetscCheck(flg, PetscObjectComm((PetscObject)smooth), PETSC_ERR_ORDER, "Pmat not set"); 23732fe681dSStefano Zampini PetscCall(KSPGetOperators(smooth, NULL, &A)); 23832fe681dSStefano Zampini } 23932fe681dSStefano Zampini 24032fe681dSStefano Zampini /* Setup (also setup smoother here) */ 24132fe681dSStefano Zampini if (!pc->setupcalled) PetscCall(KSPSetFromOptions(smooth)); 24232fe681dSStefano Zampini PetscCall(KSPSetUp(smooth)); 24332fe681dSStefano Zampini PetscCall(KSPSetUpOnBlocks(smooth)); 24432fe681dSStefano Zampini PetscCall(PCMGGDSWSetUp(pc, l, dm, smooth, Nc, A, &ns, &sA_IG, &sksp, &sI, &sG, &sGf, &sGi, &sGiM)); 24532fe681dSStefano Zampini 24632fe681dSStefano Zampini /* Number GDSW basis functions */ 24732fe681dSStefano Zampini PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A), ns, sGi, &Gidx)); 24832fe681dSStefano Zampini PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A), ns, sGiM, &GidxMult)); 24932fe681dSStefano Zampini PetscCall(ISRenumber(Gidx, GidxMult, &NG, &cG)); 25032fe681dSStefano Zampini PetscCall(ISDestroy(&Gidx)); 25132fe681dSStefano Zampini 25232fe681dSStefano Zampini /* Detect column block size */ 25332fe681dSStefano Zampini PetscCall(ISGetMinMax(GidxMult, &cbs[0], &cbs[1])); 25432fe681dSStefano Zampini PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)A), cbs, cbs)); 25532fe681dSStefano Zampini PetscCall(ISDestroy(&GidxMult)); 25632fe681dSStefano Zampini 25732fe681dSStefano Zampini /* Construct global interpolation matrix */ 25832fe681dSStefano Zampini PetscCall(MatGetLocalSize(A, NULL, &n)); 25932fe681dSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator)); 26032fe681dSStefano Zampini PetscCall(MatSetSizes(preallocator, n, PETSC_DECIDE, PETSC_DECIDE, NG)); 26132fe681dSStefano Zampini PetscCall(MatSetType(preallocator, MATPREALLOCATOR)); 26232fe681dSStefano Zampini PetscCall(MatSetUp(preallocator)); 26332fe681dSStefano Zampini PetscCall(ISGetIndices(cG, &cidx)); 26432fe681dSStefano Zampini for (i = 0, c = 0; i < ns; i++) { 26532fe681dSStefano Zampini const PetscInt *ri, *rg; 26632fe681dSStefano Zampini PetscInt nri, nrg, ncg; 26732fe681dSStefano Zampini 26832fe681dSStefano Zampini PetscCall(ISGetLocalSize(sI[i], &nri)); 26932fe681dSStefano Zampini PetscCall(ISGetLocalSize(sG[i], &nrg)); 27032fe681dSStefano Zampini PetscCall(ISGetIndices(sI[i], &ri)); 27132fe681dSStefano Zampini PetscCall(ISGetIndices(sG[i], &rg)); 27232fe681dSStefano Zampini PetscCall(MatGetSize(sGf[i], NULL, &ncg)); 27332fe681dSStefano Zampini PetscCall(MatSetValues(preallocator, nri, ri, ncg, cidx + c, NULL, INSERT_VALUES)); 27432fe681dSStefano Zampini PetscCall(MatSetValues(preallocator, nrg, rg, ncg, cidx + c, NULL, INSERT_VALUES)); 27532fe681dSStefano Zampini PetscCall(ISRestoreIndices(sI[i], &ri)); 27632fe681dSStefano Zampini PetscCall(ISRestoreIndices(sG[i], &rg)); 27732fe681dSStefano Zampini } 27832fe681dSStefano Zampini PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY)); 27932fe681dSStefano Zampini PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY)); 28032fe681dSStefano Zampini 28132fe681dSStefano Zampini ptype = MATAIJ; 28232fe681dSStefano Zampini if (PetscDefined(HAVE_DEVICE)) { 28332fe681dSStefano Zampini PetscCall(MatBoundToCPU(A, &flg)); 28432fe681dSStefano Zampini if (!flg) { 28532fe681dSStefano Zampini VecType vtype; 286bbcf679cSJacob Faibussowitsch char *found = NULL; 28732fe681dSStefano Zampini 28832fe681dSStefano Zampini PetscCall(MatGetVecType(A, &vtype)); 28932fe681dSStefano Zampini PetscCall(PetscStrstr(vtype, "cuda", &found)); 29032fe681dSStefano Zampini if (found) ptype = MATAIJCUSPARSE; 29132fe681dSStefano Zampini } 29232fe681dSStefano Zampini } 29332fe681dSStefano Zampini PetscCall(MatCreate(PetscObjectComm((PetscObject)A), cspace)); 29432fe681dSStefano Zampini PetscCall(MatSetSizes(*cspace, n, PETSC_DECIDE, PETSC_DECIDE, NG)); 29532fe681dSStefano Zampini PetscCall(MatSetType(*cspace, ptype)); 29632fe681dSStefano Zampini PetscCall(MatGetBlockSizes(A, NULL, &rbs)); 29732fe681dSStefano Zampini PetscCall(MatSetBlockSizes(*cspace, rbs, cbs[0] == cbs[1] ? cbs[0] : 1)); 29832fe681dSStefano Zampini PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_FALSE, *cspace)); 29932fe681dSStefano Zampini PetscCall(MatDestroy(&preallocator)); 30032fe681dSStefano Zampini PetscCall(MatSetOption(*cspace, MAT_ROW_ORIENTED, PETSC_FALSE)); 30132fe681dSStefano Zampini 30232fe681dSStefano Zampini for (i = 0, c = 0; i < ns; i++) { 30332fe681dSStefano Zampini Mat X, Y; 30432fe681dSStefano Zampini const PetscScalar *v; 30532fe681dSStefano Zampini const PetscInt *ri, *rg; 30632fe681dSStefano Zampini PetscInt nri, nrg, ncg; 30732fe681dSStefano Zampini 308fb842aefSJose E. Roman PetscCall(MatMatMult(sA_IG[i], sGf[i], MAT_INITIAL_MATRIX, PETSC_CURRENT, &Y)); 30932fe681dSStefano Zampini PetscCall(MatScale(Y, -1.0)); 31032fe681dSStefano Zampini PetscCall(MatDuplicate(Y, MAT_DO_NOT_COPY_VALUES, &X)); 31132fe681dSStefano Zampini PetscCall(KSPMatSolve(sksp[i], Y, X)); 31232fe681dSStefano Zampini 31332fe681dSStefano Zampini PetscCall(ISGetLocalSize(sI[i], &nri)); 31432fe681dSStefano Zampini PetscCall(ISGetLocalSize(sG[i], &nrg)); 31532fe681dSStefano Zampini PetscCall(ISGetIndices(sI[i], &ri)); 31632fe681dSStefano Zampini PetscCall(ISGetIndices(sG[i], &rg)); 31732fe681dSStefano Zampini PetscCall(MatGetSize(sGf[i], NULL, &ncg)); 31832fe681dSStefano Zampini 31932fe681dSStefano Zampini PetscCall(MatDenseGetArrayRead(X, &v)); 32032fe681dSStefano Zampini PetscCall(MatSetValues(*cspace, nri, ri, ncg, cidx + c, v, INSERT_VALUES)); 32132fe681dSStefano Zampini PetscCall(MatDenseRestoreArrayRead(X, &v)); 32232fe681dSStefano Zampini PetscCall(MatDenseGetArrayRead(sGf[i], &v)); 32332fe681dSStefano Zampini PetscCall(MatSetValues(*cspace, nrg, rg, ncg, cidx + c, v, INSERT_VALUES)); 32432fe681dSStefano Zampini PetscCall(MatDenseRestoreArrayRead(sGf[i], &v)); 32532fe681dSStefano Zampini PetscCall(ISRestoreIndices(sI[i], &ri)); 32632fe681dSStefano Zampini PetscCall(ISRestoreIndices(sG[i], &rg)); 32732fe681dSStefano Zampini PetscCall(MatDestroy(&Y)); 32832fe681dSStefano Zampini PetscCall(MatDestroy(&X)); 32932fe681dSStefano Zampini } 33032fe681dSStefano Zampini PetscCall(ISRestoreIndices(cG, &cidx)); 33132fe681dSStefano Zampini PetscCall(ISDestroy(&cG)); 33232fe681dSStefano Zampini PetscCall(MatAssemblyBegin(*cspace, MAT_FINAL_ASSEMBLY)); 33332fe681dSStefano Zampini PetscCall(MatAssemblyEnd(*cspace, MAT_FINAL_ASSEMBLY)); 33432fe681dSStefano Zampini 33532fe681dSStefano Zampini for (i = 0; i < ns; i++) { 33632fe681dSStefano Zampini PetscCall(KSPDestroy(&sksp[i])); 33732fe681dSStefano Zampini PetscCall(ISDestroy(&sI[i])); 33832fe681dSStefano Zampini PetscCall(ISDestroy(&sG[i])); 33932fe681dSStefano Zampini PetscCall(ISDestroy(&sGi[i])); 34032fe681dSStefano Zampini PetscCall(ISDestroy(&sGiM[i])); 34132fe681dSStefano Zampini PetscCall(MatDestroy(&sGf[i])); 34232fe681dSStefano Zampini PetscCall(MatDestroy(&sA_IG[i])); 34332fe681dSStefano Zampini } 34432fe681dSStefano Zampini PetscCall(PetscFree(sksp)); 34532fe681dSStefano Zampini PetscCall(PetscFree(sI)); 34632fe681dSStefano Zampini PetscCall(PetscFree(sG)); 34732fe681dSStefano Zampini PetscCall(PetscFree(sGi)); 34832fe681dSStefano Zampini PetscCall(PetscFree(sGiM)); 34932fe681dSStefano Zampini PetscCall(PetscFree(sGf)); 35032fe681dSStefano Zampini PetscCall(PetscFree(sA_IG)); 3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35232fe681dSStefano Zampini } 353