1 #include <petsc/private/pcmgimpl.h> 2 #include <petsc/private/pcbddcimpl.h> 3 #include <petsc/private/pcbddcprivateimpl.h> 4 5 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) { 6 KSP *sksp; 7 PC pcbddc = NULL, smoothpc; 8 PC_BDDC *ipcbddc; 9 PC_IS *ipcis; 10 Mat *sA_IG, *sGf, cmat, lA; 11 ISLocalToGlobalMapping l2g; 12 IS *sI, *sG, *sGi, *sGiM, cref; 13 PCBDDCSubSchurs sub_schurs = NULL; 14 PCBDDCGraph graph; 15 const char *prefix; 16 const PetscScalar *tdata; 17 PetscScalar *data, *cdata; 18 PetscReal tol = 0.0, otol; 19 const PetscInt *ia, *ja; 20 PetscInt *ccii, *cridx; 21 PetscInt i, j, ngct, ng, dbg = 0, odbg, minmax[2] = {0, PETSC_MAX_INT}, ominmax[2], vsize; 22 PetscBool flg, userdefined = PETSC_TRUE, reuse_solver = PETSC_TRUE, reduced = PETSC_FALSE; 23 24 PetscFunctionBegin; 25 PetscCall(MatGetBlockSize(A, &vsize)); 26 PetscCall(KSPGetOptionsPrefix(smooth, &prefix)); 27 PetscOptionsBegin(PetscObjectComm((PetscObject)smooth), prefix, "GDSW options", "PC"); 28 PetscCall(PetscOptionsReal("-gdsw_tolerance", "Tolerance for eigenvalue problem", NULL, tol, &tol, NULL)); 29 PetscCall(PetscOptionsBool("-gdsw_userdefined", "Use user-defined functions in addition to those adaptively generated", NULL, userdefined, &userdefined, NULL)); 30 PetscCall(PetscOptionsIntArray("-gdsw_minmax", "Minimum and maximum number of basis functions per connected component for adaptive GDSW", NULL, minmax, (i = 2, &i), NULL)); 31 PetscCall(PetscOptionsInt("-gdsw_vertex_size", "Connected components smaller or equal to vertex size will be considered as vertices", NULL, vsize, &vsize, NULL)); 32 PetscCall(PetscOptionsBool("-gdsw_reuse", "Reuse interior solver from Schur complement computations", NULL, reuse_solver, &reuse_solver, NULL)); 33 PetscCall(PetscOptionsBool("-gdsw_reduced", "Reduced GDSW", NULL, reduced, &reduced, NULL)); 34 PetscCall(PetscOptionsInt("-gdsw_debug", "Debug output", NULL, dbg, &dbg, NULL)); 35 PetscOptionsEnd(); 36 37 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATIS, &flg)); 38 if (!flg) { 39 MatNullSpace nnsp; 40 41 PetscCall(MatGetNearNullSpace(A, &nnsp)); 42 PetscObjectReference((PetscObject)nnsp); 43 PetscCall(MatConvert(A, MATIS, MAT_INITIAL_MATRIX, &A)); 44 PetscCall(MatSetNearNullSpace(A, nnsp)); 45 PetscCall(MatNullSpaceDestroy(&nnsp)); 46 } else PetscCall(PetscObjectReference((PetscObject)A)); 47 48 /* TODO Multi sub */ 49 *ns = 1; 50 PetscCall(PetscMalloc1(*ns, &sA_IG)); 51 PetscCall(PetscMalloc1(*ns, &sksp)); 52 PetscCall(PetscMalloc1(*ns, &sI)); 53 PetscCall(PetscMalloc1(*ns, &sG)); 54 PetscCall(PetscMalloc1(*ns, &sGf)); 55 PetscCall(PetscMalloc1(*ns, &sGi)); 56 PetscCall(PetscMalloc1(*ns, &sGiM)); 57 *sA_IG_n = sA_IG; 58 *sksp_n = sksp; 59 *sI_n = sI; 60 *sG_n = sG; 61 *sGf_n = sGf; 62 *sGi_n = sGi; 63 *sGiM_n = sGiM; 64 65 /* submatrices and solvers */ 66 PetscCall(KSPGetPC(smooth, &smoothpc)); 67 PetscCall(PetscObjectTypeCompareAny((PetscObject)smoothpc, &flg, PCBDDC, "")); 68 if (!flg) { 69 Mat smoothA; 70 71 PetscCall(PCGetOperators(smoothpc, &smoothA, NULL)); 72 PetscCall(PCCreate(PetscObjectComm((PetscObject)A), &pcbddc)); 73 PetscCall(PCSetType(pcbddc, PCBDDC)); 74 PetscCall(PCSetOperators(pcbddc, smoothA, A)); 75 PetscCall(PCISSetUp(pcbddc, PETSC_TRUE, PETSC_FALSE)); 76 } else { 77 PetscCall(PetscObjectReference((PetscObject)smoothpc)); 78 pcbddc = smoothpc; 79 } 80 ipcis = (PC_IS *)pcbddc->data; 81 ipcbddc = (PC_BDDC *)pcbddc->data; 82 PetscCall(PetscObjectReference((PetscObject)ipcis->A_IB)); 83 PetscCall(PetscObjectReference((PetscObject)ipcis->is_I_global)); 84 PetscCall(PetscObjectReference((PetscObject)ipcis->is_B_global)); 85 sA_IG[0] = ipcis->A_IB; 86 sI[0] = ipcis->is_I_global; 87 sG[0] = ipcis->is_B_global; 88 89 PetscCall(KSPCreate(PetscObjectComm((PetscObject)ipcis->A_II), &sksp[0])); 90 PetscCall(KSPSetOperators(sksp[0], ipcis->A_II, ipcis->pA_II)); 91 PetscCall(KSPSetOptionsPrefix(sksp[0], prefix)); 92 PetscCall(KSPAppendOptionsPrefix(sksp[0], "gdsw_")); 93 PetscCall(KSPSetFromOptions(sksp[0])); 94 95 /* analyze interface */ 96 PetscCall(MatISGetLocalMat(A, &lA)); 97 graph = ipcbddc->mat_graph; 98 if (!flg) { 99 PetscInt N; 100 101 PetscCall(MatISGetLocalToGlobalMapping(A, &l2g, NULL)); 102 PetscCall(MatGetSize(A, &N, NULL)); 103 graph->commsizelimit = 0; /* don't use the COMM_SELF variant of the graph */ 104 PetscCall(PCBDDCGraphInit(graph, l2g, N, PETSC_MAX_INT)); 105 PetscCall(MatGetRowIJ(lA, 0, PETSC_TRUE, PETSC_FALSE, &graph->nvtxs_csr, (const PetscInt **)&graph->xadj, (const PetscInt **)&graph->adjncy, &flg)); 106 PetscCall(PCBDDCGraphSetUp(graph, vsize, NULL, NULL, 0, NULL, NULL)); 107 PetscCall(MatRestoreRowIJ(lA, 0, PETSC_TRUE, PETSC_FALSE, &graph->nvtxs_csr, (const PetscInt **)&graph->xadj, (const PetscInt **)&graph->adjncy, &flg)); 108 PetscCall(PCBDDCGraphComputeConnectedComponents(graph)); 109 } 110 l2g = graph->l2gmap; 111 if (reduced) { 112 PetscContainer gcand; 113 PCBDDCGraphCandidates cand; 114 PetscErrorCode (*rgdsw)(DM, PetscInt *, IS **); 115 116 PetscCall(PetscObjectQueryFunction((PetscObject)dm, "DMComputeLocalRGDSWSets", &rgdsw)); 117 PetscCheck(rgdsw, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not supported"); 118 PetscCall(PetscNew(&cand)); 119 PetscCall((*rgdsw)(dm, &cand->nfc, &cand->Faces)); 120 /* filter interior (if any) and guarantee IS are ordered by global numbering */ 121 for (i = 0; i < cand->nfc; i++) { 122 IS is, is2; 123 124 PetscCall(ISLocalToGlobalMappingApplyIS(l2g, cand->Faces[i], &is)); 125 PetscCall(ISDestroy(&cand->Faces[i])); 126 PetscCall(ISSort(is)); 127 PetscCall(ISGlobalToLocalMappingApplyIS(l2g, IS_GTOLM_DROP, is, &is2)); 128 PetscCall(ISDestroy(&is)); 129 PetscCall(ISGlobalToLocalMappingApplyIS(ipcis->BtoNmap, IS_GTOLM_DROP, is2, &is)); 130 PetscCall(ISDestroy(&is2)); 131 PetscCall(ISLocalToGlobalMappingApplyIS(ipcis->BtoNmap, is, &cand->Faces[i])); 132 PetscCall(ISDestroy(&is)); 133 } 134 PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &gcand)); 135 PetscCall(PetscContainerSetPointer(gcand, cand)); 136 PetscCall(PetscContainerSetUserDestroy(gcand, PCBDDCDestroyGraphCandidatesIS)); 137 PetscCall(PetscObjectCompose((PetscObject)l2g, "_PCBDDCGraphCandidatesIS", (PetscObject)gcand)); 138 PetscCall(PetscContainerDestroy(&gcand)); 139 } 140 141 /* interface functions */ 142 otol = ipcbddc->adaptive_threshold[1]; 143 odbg = ipcbddc->dbg_flag; 144 ominmax[0] = ipcbddc->adaptive_nmin; 145 ominmax[1] = ipcbddc->adaptive_nmax; 146 ipcbddc->adaptive_threshold[1] = tol; 147 ipcbddc->dbg_flag = dbg; 148 ipcbddc->adaptive_nmin = minmax[0]; 149 ipcbddc->adaptive_nmax = minmax[1]; 150 if (tol != 0.0) { /* adaptive */ 151 Mat lS; 152 153 PetscCall(MatCreateSchurComplement(ipcis->A_II, ipcis->pA_II, ipcis->A_IB, ipcis->A_BI, ipcis->A_BB, &lS)); 154 PetscCall(KSPGetOptionsPrefix(sksp[0], &prefix)); 155 PetscCall(PCBDDCSubSchursCreate(&sub_schurs)); 156 PetscCall(PCBDDCSubSchursInit(sub_schurs, prefix, ipcis->is_I_local, ipcis->is_B_local, graph, ipcis->BtoNmap, PETSC_FALSE, PETSC_TRUE)); 157 if (userdefined) PetscCall(PCBDDCComputeFakeChange(pcbddc, PETSC_FALSE, graph, NULL, &cmat, &cref, NULL, &flg)); 158 else { 159 cmat = NULL; 160 cref = NULL; 161 } 162 PetscCall(PCBDDCSubSchursSetUp(sub_schurs, lA, lS, PETSC_TRUE, NULL, NULL, -1, NULL, PETSC_TRUE, reuse_solver, PETSC_FALSE, 0, NULL, NULL, cmat, cref)); 163 PetscCall(MatDestroy(&lS)); 164 PetscCall(MatDestroy(&cmat)); 165 PetscCall(ISDestroy(&cref)); 166 if (sub_schurs->reuse_solver) { 167 PetscCall(KSPSetPC(sksp[0], sub_schurs->reuse_solver->interior_solver)); 168 PetscCall(PCDestroy(&sub_schurs->reuse_solver->interior_solver)); 169 sub_schurs->reuse_solver = NULL; 170 } 171 } 172 PetscCall(PCBDDCComputeFakeChange(pcbddc, PETSC_TRUE, graph, sub_schurs, &cmat, &cref, &sGiM[0], NULL)); 173 PetscCall(PCBDDCSubSchursDestroy(&sub_schurs)); 174 ipcbddc->adaptive_threshold[1] = otol; 175 ipcbddc->dbg_flag = odbg; 176 ipcbddc->adaptive_nmin = ominmax[0]; 177 ipcbddc->adaptive_nmax = ominmax[1]; 178 179 PetscCall(ISLocalToGlobalMappingApplyIS(l2g, cref, &sGi[0])); 180 PetscCall(ISDestroy(&cref)); 181 182 PetscCall(MatSeqAIJGetArrayRead(cmat, &tdata)); 183 PetscCall(MatGetRowIJ(cmat, 0, PETSC_FALSE, PETSC_FALSE, &ngct, &ia, &ja, &flg)); 184 PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in MatGetRowIJ"); 185 186 PetscCall(PetscMalloc1(ngct + 1, &ccii)); 187 PetscCall(PetscMalloc1(ia[ngct], &cridx)); 188 PetscCall(PetscMalloc1(ia[ngct], &cdata)); 189 190 PetscCall(PetscArraycpy(ccii, ia, ngct + 1)); 191 PetscCall(PetscArraycpy(cdata, tdata, ia[ngct])); 192 PetscCall(ISGlobalToLocalMappingApply(ipcis->BtoNmap, IS_GTOLM_DROP, ia[ngct], ja, &i, cridx)); 193 PetscCheck(i == ia[ngct], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in G2L"); 194 195 PetscCall(MatRestoreRowIJ(cmat, 0, PETSC_FALSE, PETSC_FALSE, &i, &ia, &ja, &flg)); 196 PetscCall(MatSeqAIJRestoreArrayRead(cmat, &tdata)); 197 PetscCall(MatDestroy(&cmat)); 198 199 /* populate dense matrix */ 200 PetscCall(ISGetLocalSize(sG[0], &ng)); 201 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, ng, ngct, NULL, &sGf[0])); 202 PetscCall(MatDenseGetArrayWrite(sGf[0], &data)); 203 for (i = 0; i < ngct; i++) 204 for (j = ccii[i]; j < ccii[i + 1]; j++) data[ng * i + cridx[j]] = cdata[j]; 205 PetscCall(MatDenseRestoreArrayWrite(sGf[0], &data)); 206 207 PetscCall(PetscFree(cdata)); 208 PetscCall(PetscFree(ccii)); 209 PetscCall(PetscFree(cridx)); 210 PetscCall(PCDestroy(&pcbddc)); 211 PetscCall(MatDestroy(&A)); 212 PetscFunctionReturn(0); 213 } 214 215 PetscErrorCode PCMGGDSWCreateCoarseSpace_Private(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat guess, Mat *cspace) { 216 KSP *sksp; 217 Mat A, *sA_IG, *sGf, preallocator; 218 IS Gidx, GidxMult, cG; 219 IS *sI, *sG, *sGi, *sGiM; 220 const PetscInt *cidx; 221 PetscInt NG, ns, n, i, c, rbs, cbs[2]; 222 PetscBool flg; 223 MatType ptype; 224 225 PetscFunctionBegin; 226 *cspace = NULL; 227 if (!l) PetscFunctionReturn(0); 228 if (pc->useAmat) { 229 PetscCall(KSPGetOperatorsSet(smooth, &flg, NULL)); 230 PetscCheck(flg, PetscObjectComm((PetscObject)smooth), PETSC_ERR_ORDER, "Amat not set"); 231 PetscCall(KSPGetOperators(smooth, &A, NULL)); 232 } else { 233 PetscCall(KSPGetOperatorsSet(smooth, NULL, &flg)); 234 PetscCheck(flg, PetscObjectComm((PetscObject)smooth), PETSC_ERR_ORDER, "Pmat not set"); 235 PetscCall(KSPGetOperators(smooth, NULL, &A)); 236 } 237 238 /* Setup (also setup smoother here) */ 239 if (!pc->setupcalled) PetscCall(KSPSetFromOptions(smooth)); 240 PetscCall(KSPSetUp(smooth)); 241 PetscCall(KSPSetUpOnBlocks(smooth)); 242 PetscCall(PCMGGDSWSetUp(pc, l, dm, smooth, Nc, A, &ns, &sA_IG, &sksp, &sI, &sG, &sGf, &sGi, &sGiM)); 243 244 /* Number GDSW basis functions */ 245 PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A), ns, sGi, &Gidx)); 246 PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A), ns, sGiM, &GidxMult)); 247 PetscCall(ISRenumber(Gidx, GidxMult, &NG, &cG)); 248 PetscCall(ISDestroy(&Gidx)); 249 250 /* Detect column block size */ 251 PetscCall(ISGetMinMax(GidxMult, &cbs[0], &cbs[1])); 252 PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)A), cbs, cbs)); 253 PetscCall(ISDestroy(&GidxMult)); 254 255 /* Construct global interpolation matrix */ 256 PetscCall(MatGetLocalSize(A, NULL, &n)); 257 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator)); 258 PetscCall(MatSetSizes(preallocator, n, PETSC_DECIDE, PETSC_DECIDE, NG)); 259 PetscCall(MatSetType(preallocator, MATPREALLOCATOR)); 260 PetscCall(MatSetUp(preallocator)); 261 PetscCall(ISGetIndices(cG, &cidx)); 262 for (i = 0, c = 0; i < ns; i++) { 263 const PetscInt *ri, *rg; 264 PetscInt nri, nrg, ncg; 265 266 PetscCall(ISGetLocalSize(sI[i], &nri)); 267 PetscCall(ISGetLocalSize(sG[i], &nrg)); 268 PetscCall(ISGetIndices(sI[i], &ri)); 269 PetscCall(ISGetIndices(sG[i], &rg)); 270 PetscCall(MatGetSize(sGf[i], NULL, &ncg)); 271 PetscCall(MatSetValues(preallocator, nri, ri, ncg, cidx + c, NULL, INSERT_VALUES)); 272 PetscCall(MatSetValues(preallocator, nrg, rg, ncg, cidx + c, NULL, INSERT_VALUES)); 273 PetscCall(ISRestoreIndices(sI[i], &ri)); 274 PetscCall(ISRestoreIndices(sG[i], &rg)); 275 } 276 PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY)); 277 PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY)); 278 279 ptype = MATAIJ; 280 if (PetscDefined(HAVE_DEVICE)) { 281 PetscCall(MatBoundToCPU(A, &flg)); 282 if (!flg) { 283 VecType vtype; 284 char *found; 285 286 PetscCall(MatGetVecType(A, &vtype)); 287 PetscCall(PetscStrstr(vtype, "cuda", &found)); 288 if (found) ptype = MATAIJCUSPARSE; 289 } 290 } 291 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), cspace)); 292 PetscCall(MatSetSizes(*cspace, n, PETSC_DECIDE, PETSC_DECIDE, NG)); 293 PetscCall(MatSetType(*cspace, ptype)); 294 PetscCall(MatGetBlockSizes(A, NULL, &rbs)); 295 PetscCall(MatSetBlockSizes(*cspace, rbs, cbs[0] == cbs[1] ? cbs[0] : 1)); 296 PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_FALSE, *cspace)); 297 PetscCall(MatDestroy(&preallocator)); 298 PetscCall(MatSetOption(*cspace, MAT_ROW_ORIENTED, PETSC_FALSE)); 299 300 for (i = 0, c = 0; i < ns; i++) { 301 Mat X, Y; 302 const PetscScalar *v; 303 const PetscInt *ri, *rg; 304 PetscInt nri, nrg, ncg; 305 306 PetscCall(MatMatMult(sA_IG[i], sGf[i], MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Y)); 307 PetscCall(MatScale(Y, -1.0)); 308 PetscCall(MatDuplicate(Y, MAT_DO_NOT_COPY_VALUES, &X)); 309 PetscCall(KSPMatSolve(sksp[i], Y, X)); 310 311 PetscCall(ISGetLocalSize(sI[i], &nri)); 312 PetscCall(ISGetLocalSize(sG[i], &nrg)); 313 PetscCall(ISGetIndices(sI[i], &ri)); 314 PetscCall(ISGetIndices(sG[i], &rg)); 315 PetscCall(MatGetSize(sGf[i], NULL, &ncg)); 316 317 PetscCall(MatDenseGetArrayRead(X, &v)); 318 PetscCall(MatSetValues(*cspace, nri, ri, ncg, cidx + c, v, INSERT_VALUES)); 319 PetscCall(MatDenseRestoreArrayRead(X, &v)); 320 PetscCall(MatDenseGetArrayRead(sGf[i], &v)); 321 PetscCall(MatSetValues(*cspace, nrg, rg, ncg, cidx + c, v, INSERT_VALUES)); 322 PetscCall(MatDenseRestoreArrayRead(sGf[i], &v)); 323 PetscCall(ISRestoreIndices(sI[i], &ri)); 324 PetscCall(ISRestoreIndices(sG[i], &rg)); 325 PetscCall(MatDestroy(&Y)); 326 PetscCall(MatDestroy(&X)); 327 } 328 PetscCall(ISRestoreIndices(cG, &cidx)); 329 PetscCall(ISDestroy(&cG)); 330 PetscCall(MatAssemblyBegin(*cspace, MAT_FINAL_ASSEMBLY)); 331 PetscCall(MatAssemblyEnd(*cspace, MAT_FINAL_ASSEMBLY)); 332 333 for (i = 0; i < ns; i++) { 334 PetscCall(KSPDestroy(&sksp[i])); 335 PetscCall(ISDestroy(&sI[i])); 336 PetscCall(ISDestroy(&sG[i])); 337 PetscCall(ISDestroy(&sGi[i])); 338 PetscCall(ISDestroy(&sGiM[i])); 339 PetscCall(MatDestroy(&sGf[i])); 340 PetscCall(MatDestroy(&sA_IG[i])); 341 } 342 PetscCall(PetscFree(sksp)); 343 PetscCall(PetscFree(sI)); 344 PetscCall(PetscFree(sG)); 345 PetscCall(PetscFree(sGi)); 346 PetscCall(PetscFree(sGiM)); 347 PetscCall(PetscFree(sGf)); 348 PetscCall(PetscFree(sA_IG)); 349 PetscFunctionReturn(0); 350 } 351