1 #include <petsc/private/pcpatchimpl.h> /*I "petscpc.h" I*/ 2 #include <petscdmplex.h> 3 #include <petscsf.h> 4 #include <petscbt.h> 5 6 PetscLogEvent PC_Patch_CreatePatches, PC_Patch_ComputeOp, PC_Patch_Solve, PC_Patch_Scatter, PC_Patch_Apply, PC_Patch_Prealloc; 7 8 static PetscErrorCode PCPatchConstruct_Star(void *vpatch, DM dm, PetscInt point, PetscHashI ht) 9 { 10 PetscInt starSize; 11 PetscInt *star = NULL, si; 12 PetscErrorCode ierr; 13 14 PetscFunctionBegin; 15 PetscHashIClear(ht); 16 /* To start with, add the point we care about */ 17 PetscHashIAdd(ht, point, 0); 18 /* Loop over all the points that this point connects to */ 19 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 20 for (si = 0; si < starSize; si += 2) {PetscHashIAdd(ht, star[si], 0);} 21 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 22 PetscFunctionReturn(0); 23 } 24 25 static PetscErrorCode PCPatchConstruct_Vanka(void *vpatch, DM dm, PetscInt point, PetscHashI ht) 26 { 27 PC_PATCH *patch = (PC_PATCH *) vpatch; 28 PetscInt starSize; 29 PetscInt *star = NULL; 30 PetscBool shouldIgnore = PETSC_FALSE; 31 PetscInt cStart, cEnd, iStart, iEnd, si; 32 PetscErrorCode ierr; 33 34 PetscFunctionBegin; 35 PetscHashIClear(ht); 36 /* To start with, add the point we care about */ 37 PetscHashIAdd(ht, point, 0); 38 /* Should we ignore any points of a certain dimension? */ 39 if (patch->vankadim >= 0) { 40 shouldIgnore = PETSC_TRUE; 41 ierr = DMPlexGetDepthStratum(dm, patch->vankadim, &iStart, &iEnd);CHKERRQ(ierr); 42 } 43 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 44 /* Loop over all the cells that this point connects to */ 45 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 46 for (si = 0; si < starSize; si += 2) { 47 const PetscInt cell = star[si]; 48 PetscInt closureSize; 49 PetscInt *closure = NULL, ci; 50 51 if (cell < cStart || cell >= cEnd) continue; 52 /* now loop over all entities in the closure of that cell */ 53 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 54 for (ci = 0; ci < closureSize; ci += 2) { 55 const PetscInt newpoint = closure[ci]; 56 57 /* We've been told to ignore entities of this type.*/ 58 if (shouldIgnore && newpoint >= iStart && newpoint < iEnd) continue; 59 PetscHashIAdd(ht, newpoint, 0); 60 } 61 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 62 } 63 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 64 PetscFunctionReturn(0); 65 } 66 67 /* The user's already set the patches in patch->userIS. Build the hash tables */ 68 static PetscErrorCode PCPatchConstruct_User(void *vpatch, DM dm, PetscInt point, PetscHashI ht) 69 { 70 PC_PATCH *patch = (PC_PATCH *) vpatch; 71 IS patchis = patch->userIS[point]; 72 PetscInt n; 73 const PetscInt *patchdata; 74 PetscInt pStart, pEnd, i; 75 PetscErrorCode ierr; 76 77 PetscFunctionBegin; 78 PetscHashIClear(ht); 79 ierr = DMPlexGetChart(dm, &pStart, &pEnd); 80 ierr = ISGetLocalSize(patchis, &n);CHKERRQ(ierr); 81 ierr = ISGetIndices(patchis, &patchdata);CHKERRQ(ierr); 82 for (i = 0; i < n; ++i) { 83 const PetscInt ownedpoint = patchdata[i]; 84 85 if (ownedpoint < pStart || ownedpoint >= pEnd) { 86 SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D was not in [%D, %D)", ownedpoint, pStart, pEnd); 87 } 88 PetscHashIAdd(ht, ownedpoint, 0); 89 } 90 ierr = ISRestoreIndices(patchis, &patchdata);CHKERRQ(ierr); 91 PetscFunctionReturn(0); 92 } 93 94 static PetscErrorCode PCPatchCreateDefaultSF_Private(PC pc, PetscInt n, const PetscSF *sf, const PetscInt *bs) 95 { 96 PC_PATCH *patch = (PC_PATCH *) pc->data; 97 PetscInt i; 98 PetscErrorCode ierr; 99 100 PetscFunctionBegin; 101 if (n == 1 && bs[0] == 1) { 102 patch->defaultSF = sf[0]; 103 ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr); 104 } else { 105 PetscInt allRoots = 0, allLeaves = 0; 106 PetscInt leafOffset = 0; 107 PetscInt *ilocal = NULL; 108 PetscSFNode *iremote = NULL; 109 PetscInt *remoteOffsets = NULL; 110 PetscInt index = 0; 111 PetscHashI rankToIndex; 112 PetscInt numRanks = 0; 113 PetscSFNode *remote = NULL; 114 PetscSF rankSF; 115 PetscInt *ranks = NULL; 116 PetscInt *offsets = NULL; 117 MPI_Datatype contig; 118 PetscHashI ranksUniq; 119 120 /* First figure out how many dofs there are in the concatenated numbering. 121 * allRoots: number of owned global dofs; 122 * allLeaves: number of visible dofs (global + ghosted). 123 */ 124 for (i = 0; i < n; ++i) { 125 PetscInt nroots, nleaves; 126 127 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 128 allRoots += nroots * bs[i]; 129 allLeaves += nleaves * bs[i]; 130 } 131 ierr = PetscMalloc1(allLeaves, &ilocal);CHKERRQ(ierr); 132 ierr = PetscMalloc1(allLeaves, &iremote);CHKERRQ(ierr); 133 /* Now build an SF that just contains process connectivity. */ 134 PetscHashICreate(ranksUniq); 135 for (i = 0; i < n; ++i) { 136 const PetscMPIInt *ranks = NULL; 137 PetscInt nranks, j; 138 139 ierr = PetscSFSetUp(sf[i]);CHKERRQ(ierr); 140 ierr = PetscSFGetRanks(sf[i], &nranks, &ranks, NULL, NULL, NULL);CHKERRQ(ierr); 141 /* These are all the ranks who communicate with me. */ 142 for (j = 0; j < nranks; ++j) { 143 PetscHashIAdd(ranksUniq, (PetscInt) ranks[j], 0); 144 } 145 } 146 PetscHashISize(ranksUniq, numRanks); 147 ierr = PetscMalloc1(numRanks, &remote);CHKERRQ(ierr); 148 ierr = PetscMalloc1(numRanks, &ranks);CHKERRQ(ierr); 149 PetscHashIGetKeys(ranksUniq, &index, ranks); 150 151 PetscHashICreate(rankToIndex); 152 for (i = 0; i < numRanks; ++i) { 153 remote[i].rank = ranks[i]; 154 remote[i].index = 0; 155 PetscHashIAdd(rankToIndex, ranks[i], i); 156 } 157 ierr = PetscFree(ranks);CHKERRQ(ierr); 158 PetscHashIDestroy(ranksUniq); 159 ierr = PetscSFCreate(PetscObjectComm((PetscObject) pc), &rankSF);CHKERRQ(ierr); 160 ierr = PetscSFSetGraph(rankSF, 1, numRanks, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 161 ierr = PetscSFSetUp(rankSF);CHKERRQ(ierr); 162 /* OK, use it to communicate the root offset on the remote 163 * processes for each subspace. */ 164 ierr = PetscMalloc1(n, &offsets);CHKERRQ(ierr); 165 ierr = PetscMalloc1(n*numRanks, &remoteOffsets);CHKERRQ(ierr); 166 167 offsets[0] = 0; 168 for (i = 1; i < n; ++i) { 169 PetscInt nroots; 170 171 ierr = PetscSFGetGraph(sf[i-1], &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 172 offsets[i] = offsets[i-1] + nroots*bs[i-1]; 173 } 174 /* Offsets are the offsets on the current process of the 175 * global dof numbering for the subspaces. */ 176 ierr = MPI_Type_contiguous(n, MPIU_INT, &contig);CHKERRQ(ierr); 177 ierr = MPI_Type_commit(&contig);CHKERRQ(ierr); 178 179 ierr = PetscSFBcastBegin(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 180 ierr = PetscSFBcastEnd(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 181 ierr = MPI_Type_free(&contig);CHKERRQ(ierr); 182 ierr = PetscFree(offsets);CHKERRQ(ierr); 183 ierr = PetscSFDestroy(&rankSF);CHKERRQ(ierr); 184 /* Now remoteOffsets contains the offsets on the remote 185 * processes who communicate with me. So now we can 186 * concatenate the list of SFs into a single one. */ 187 index = 0; 188 for (i = 0; i < n; ++i) { 189 const PetscSFNode *remote = NULL; 190 const PetscInt *local = NULL; 191 PetscInt nroots, nleaves, j; 192 193 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, &local, &remote);CHKERRQ(ierr); 194 for (j = 0; j < nleaves; ++j) { 195 PetscInt rank = remote[j].rank; 196 PetscInt idx, rootOffset, k; 197 198 PetscHashIMap(rankToIndex, rank, idx); 199 if (idx == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Didn't find rank, huh?"); 200 /* Offset on given rank for ith subspace */ 201 rootOffset = remoteOffsets[n*idx + i]; 202 for (k = 0; k < bs[i]; ++k) { 203 ilocal[index] = local[j]*bs[i] + k + leafOffset; 204 iremote[index].rank = remote[j].rank; 205 iremote[index].index = remote[j].index*bs[i] + k + rootOffset; 206 ++index; 207 } 208 } 209 leafOffset += nleaves * bs[i]; 210 } 211 PetscHashIDestroy(rankToIndex); 212 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 213 ierr = PetscSFCreate(PetscObjectComm((PetscObject)pc), &patch->defaultSF);CHKERRQ(ierr); 214 ierr = PetscSFSetGraph(patch->defaultSF, allRoots, allLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 215 } 216 PetscFunctionReturn(0); 217 } 218 219 /* TODO: Docs */ 220 PetscErrorCode PCPatchSetSaveOperators(PC pc, PetscBool flg) 221 { 222 PC_PATCH *patch = (PC_PATCH *) pc->data; 223 PetscFunctionBegin; 224 patch->save_operators = flg; 225 PetscFunctionReturn(0); 226 } 227 228 /* TODO: Docs */ 229 PetscErrorCode PCPatchGetSaveOperators(PC pc, PetscBool *flg) 230 { 231 PC_PATCH *patch = (PC_PATCH *) pc->data; 232 PetscFunctionBegin; 233 *flg = patch->save_operators; 234 PetscFunctionReturn(0); 235 } 236 237 /* TODO: Docs */ 238 PetscErrorCode PCPatchSetPartitionOfUnity(PC pc, PetscBool flg) 239 { 240 PC_PATCH *patch = (PC_PATCH *) pc->data; 241 PetscFunctionBegin; 242 patch->partition_of_unity = flg; 243 PetscFunctionReturn(0); 244 } 245 246 /* TODO: Docs */ 247 PetscErrorCode PCPatchGetPartitionOfUnity(PC pc, PetscBool *flg) 248 { 249 PC_PATCH *patch = (PC_PATCH *) pc->data; 250 PetscFunctionBegin; 251 *flg = patch->partition_of_unity; 252 PetscFunctionReturn(0); 253 } 254 255 /* TODO: Docs */ 256 PetscErrorCode PCPatchSetMultiplicative(PC pc, PetscBool flg) 257 { 258 PC_PATCH *patch = (PC_PATCH *) pc->data; 259 PetscFunctionBegin; 260 patch->multiplicative = flg; 261 PetscFunctionReturn(0); 262 } 263 264 /* TODO: Docs */ 265 PetscErrorCode PCPatchGetMultiplicative(PC pc, PetscBool *flg) 266 { 267 PC_PATCH *patch = (PC_PATCH *) pc->data; 268 PetscFunctionBegin; 269 *flg = patch->multiplicative; 270 PetscFunctionReturn(0); 271 } 272 273 /* TODO: Docs */ 274 PetscErrorCode PCPatchSetSubMatType(PC pc, MatType sub_mat_type) 275 { 276 PC_PATCH *patch = (PC_PATCH *) pc->data; 277 PetscErrorCode ierr; 278 279 PetscFunctionBegin; 280 if (patch->sub_mat_type) {ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);} 281 ierr = PetscStrallocpy(sub_mat_type, (char **) &patch->sub_mat_type);CHKERRQ(ierr); 282 PetscFunctionReturn(0); 283 } 284 285 /* TODO: Docs */ 286 PetscErrorCode PCPatchGetSubMatType(PC pc, MatType *sub_mat_type) 287 { 288 PC_PATCH *patch = (PC_PATCH *) pc->data; 289 PetscFunctionBegin; 290 *sub_mat_type = patch->sub_mat_type; 291 PetscFunctionReturn(0); 292 } 293 294 /* TODO: Docs */ 295 PetscErrorCode PCPatchSetCellNumbering(PC pc, PetscSection cellNumbering) 296 { 297 PC_PATCH *patch = (PC_PATCH *) pc->data; 298 PetscErrorCode ierr; 299 300 PetscFunctionBegin; 301 patch->cellNumbering = cellNumbering; 302 ierr = PetscObjectReference((PetscObject) cellNumbering);CHKERRQ(ierr); 303 PetscFunctionReturn(0); 304 } 305 306 /* TODO: Docs */ 307 PetscErrorCode PCPatchGetCellNumbering(PC pc, PetscSection *cellNumbering) 308 { 309 PC_PATCH *patch = (PC_PATCH *) pc->data; 310 PetscFunctionBegin; 311 *cellNumbering = patch->cellNumbering; 312 PetscFunctionReturn(0); 313 } 314 315 /* TODO: Docs */ 316 PetscErrorCode PCPatchSetConstructType(PC pc, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx) 317 { 318 PC_PATCH *patch = (PC_PATCH *) pc->data; 319 320 PetscFunctionBegin; 321 patch->ctype = ctype; 322 switch (ctype) { 323 case PC_PATCH_STAR: 324 patch->patchconstructop = PCPatchConstruct_Star; 325 break; 326 case PC_PATCH_VANKA: 327 patch->patchconstructop = PCPatchConstruct_Vanka; 328 break; 329 case PC_PATCH_USER: 330 case PC_PATCH_PYTHON: 331 patch->user_patches = PETSC_TRUE; 332 patch->patchconstructop = PCPatchConstruct_User; 333 patch->userpatchconstructionop = func; 334 patch->userpatchconstructctx = ctx; 335 break; 336 default: 337 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 338 } 339 PetscFunctionReturn(0); 340 } 341 342 /* TODO: Docs */ 343 PetscErrorCode PCPatchGetConstructType(PC pc, PCPatchConstructType *ctype, PetscErrorCode (**func)(PC, PetscInt *, IS **, IS *, void *), void **ctx) 344 { 345 PC_PATCH *patch = (PC_PATCH *) pc->data; 346 347 PetscFunctionBegin; 348 *ctype = patch->ctype; 349 switch (patch->ctype) { 350 case PC_PATCH_STAR: 351 case PC_PATCH_VANKA: 352 break; 353 case PC_PATCH_USER: 354 case PC_PATCH_PYTHON: 355 *func = patch->userpatchconstructionop; 356 *ctx = patch->userpatchconstructctx; 357 break; 358 default: 359 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 360 } 361 PetscFunctionReturn(0); 362 } 363 364 /* TODO: Docs */ 365 PetscErrorCode PCPatchSetDiscretisationInfo(PC pc, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, 366 const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes) 367 { 368 PC_PATCH *patch = (PC_PATCH *) pc->data; 369 PetscSF *sfs; 370 PetscInt i; 371 PetscErrorCode ierr; 372 373 PetscFunctionBegin; 374 ierr = PetscMalloc1(nsubspaces, &sfs);CHKERRQ(ierr); 375 ierr = PetscMalloc1(nsubspaces, &patch->dofSection);CHKERRQ(ierr); 376 ierr = PetscMalloc1(nsubspaces, &patch->bs);CHKERRQ(ierr); 377 ierr = PetscMalloc1(nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr); 378 ierr = PetscMalloc1(nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr); 379 ierr = PetscMalloc1(nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr); 380 381 patch->nsubspaces = nsubspaces; 382 patch->totalDofsPerCell = 0; 383 for (i = 0; i < nsubspaces; ++i) { 384 ierr = DMGetDefaultSection(dms[i], &patch->dofSection[i]);CHKERRQ(ierr); 385 ierr = PetscObjectReference((PetscObject) patch->dofSection[i]);CHKERRQ(ierr); 386 ierr = DMGetDefaultSF(dms[i], &sfs[i]);CHKERRQ(ierr); 387 patch->bs[i] = bs[i]; 388 patch->nodesPerCell[i] = nodesPerCell[i]; 389 patch->totalDofsPerCell += nodesPerCell[i]*bs[i]; 390 patch->cellNodeMap[i] = cellNodeMap[i]; 391 patch->subspaceOffsets[i] = subspaceOffsets[i]; 392 } 393 ierr = PCPatchCreateDefaultSF_Private(pc, nsubspaces, sfs, patch->bs);CHKERRQ(ierr); 394 ierr = PetscFree(sfs);CHKERRQ(ierr); 395 396 patch->subspaceOffsets[nsubspaces] = subspaceOffsets[nsubspaces]; 397 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr); 398 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr); 399 PetscFunctionReturn(0); 400 } 401 402 /* TODO: Docs */ 403 PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, Mat, PetscInt, const PetscInt *, PetscInt, const PetscInt *, void *), void *ctx) 404 { 405 PC_PATCH *patch = (PC_PATCH *) pc->data; 406 407 PetscFunctionBegin; 408 /* User op can assume matrix is zeroed */ 409 patch->usercomputeop = func; 410 patch->usercomputectx = ctx; 411 PetscFunctionReturn(0); 412 } 413 414 /* On entry, ht contains the topological entities whose dofs we are responsible for solving for; 415 on exit, cht contains all the topological entities we need to compute their residuals. 416 In full generality this should incorporate knowledge of the sparsity pattern of the matrix; 417 here we assume a standard FE sparsity pattern.*/ 418 /* TODO: Use DMPlexGetAdjacency() */ 419 /* TODO: Look at temp buffer management for GetClosure() */ 420 static PetscErrorCode PCPatchCompleteCellPatch(DM dm, PetscHashI ht, PetscHashI cht) 421 { 422 PetscHashIIter hi; 423 PetscInt point; 424 PetscInt *star = NULL, *closure = NULL; 425 PetscInt starSize, closureSize, si, ci; 426 PetscErrorCode ierr; 427 428 PetscFunctionBegin; 429 PetscHashIClear(cht); 430 PetscHashIIterBegin(ht, hi); 431 while (!PetscHashIIterAtEnd(ht, hi)) { 432 PetscHashIIterGetKey(ht, hi, point); 433 PetscHashIIterNext(ht, hi); 434 435 /* Loop over all the cells that this point connects to */ 436 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 437 for (si = 0; si < starSize; si += 2) { 438 PetscInt ownedpoint = star[si]; 439 /* now loop over all entities in the closure of that cell */ 440 ierr = DMPlexGetTransitiveClosure(dm, ownedpoint, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 441 for (ci = 0; ci < closureSize; ci += 2) { 442 PetscInt seenpoint = closure[ci]; 443 PetscHashIAdd(cht, seenpoint, 0); 444 } 445 } 446 } 447 /* Only restore work arrays at very end. */ 448 if (closure) {ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr);} 449 if (star) {ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_FALSE, NULL, &star);CHKERRQ(ierr);} 450 PetscFunctionReturn(0); 451 } 452 453 /* Given a hash table with a set of topological entities (pts), compute the degrees of 454 freedom in global concatenated numbering on those entities. 455 For Vanka smoothing, this needs to do something special: ignore dofs of the 456 constraint subspace on entities that aren't the base entity we're building the patch 457 around. */ 458 static PetscErrorCode PCPatchGetPointDofs(PC_PATCH *patch, PetscHashI pts, PetscHashI dofs, PetscInt base, PetscInt exclude_subspace) 459 { 460 PetscHashIIter hi; 461 PetscInt ldof, loff; 462 PetscInt k, p; 463 PetscErrorCode ierr; 464 465 PetscFunctionBegin; 466 PetscHashIClear(dofs); 467 for (k = 0; k < patch->nsubspaces; ++k) { 468 PetscSection dofSection = patch->dofSection[k]; 469 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 470 PetscInt bs = patch->bs[k]; 471 PetscInt j, l; 472 473 if (k == exclude_subspace) { 474 /* only get this subspace dofs at the base entity, not any others */ 475 ierr = PetscSectionGetDof(dofSection, base, &ldof);CHKERRQ(ierr); 476 ierr = PetscSectionGetOffset(dofSection, base, &loff);CHKERRQ(ierr); 477 if (0 == ldof) continue; 478 for (j = loff; j < ldof + loff; ++j) { 479 for (l = 0; l < bs; ++l) { 480 PetscInt dof = bs*j + l + subspaceOffset; 481 PetscHashIAdd(dofs, dof, 0); 482 } 483 } 484 continue; /* skip the other dofs of this subspace */ 485 } 486 487 PetscHashIIterBegin(pts, hi); 488 while (!PetscHashIIterAtEnd(pts, hi)) { 489 PetscHashIIterGetKey(pts, hi, p); 490 PetscHashIIterNext(pts, hi); 491 ierr = PetscSectionGetDof(dofSection, p, &ldof);CHKERRQ(ierr); 492 ierr = PetscSectionGetOffset(dofSection, p, &loff);CHKERRQ(ierr); 493 if (0 == ldof) continue; 494 for (j = loff; j < ldof + loff; ++j) { 495 for (l = 0; l < bs; ++l) { 496 PetscInt dof = bs*j + l + subspaceOffset; 497 PetscHashIAdd(dofs, dof, 0); 498 } 499 } 500 } 501 } 502 PetscFunctionReturn(0); 503 } 504 505 /* Given two hash tables A and B, compute the keys in B that are not in A, and put them in C */ 506 static PetscErrorCode PCPatchComputeSetDifference_Private(PetscHashI A, PetscHashI B, PetscHashI C) 507 { 508 PetscHashIIter hi; 509 PetscInt key, val; 510 PetscBool flg; 511 512 PetscFunctionBegin; 513 PetscHashIClear(C); 514 PetscHashIIterBegin(B, hi); 515 while (!PetscHashIIterAtEnd(B, hi)) { 516 PetscHashIIterGetKeyVal(B, hi, key, val); 517 PetscHashIIterNext(B, hi); 518 PetscHashIHasKey(A, key, flg); 519 if (!flg) {PetscHashIAdd(C, key, val);} 520 } 521 PetscFunctionReturn(0); 522 } 523 524 /* 525 * PCPatchCreateCellPatches - create patches. 526 * 527 * Input Parameters: 528 * + dm - The DMPlex object defining the mesh 529 * 530 * Output Parameters: 531 * + cellCounts - Section with counts of cells around each vertex 532 * - cells - IS of the cell point indices of cells in each patch 533 */ 534 static PetscErrorCode PCPatchCreateCellPatches(PC pc) 535 { 536 PC_PATCH *patch = (PC_PATCH *) pc->data; 537 DM dm, plex; 538 DMLabel ghost; 539 PetscHashI ht, cht; 540 PetscInt *cellsArray = NULL; 541 PetscInt numCells; 542 PetscSection cellCounts; 543 PetscInt pStart, pEnd, vStart, vEnd, v, cStart, cEnd; 544 PetscBool flg; 545 PetscErrorCode ierr; 546 547 PetscFunctionBegin; 548 /* Used to keep track of the cells in the patch. */ 549 PetscHashICreate(ht); 550 PetscHashICreate(cht); 551 552 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 553 if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch PC\n"); 554 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 555 ierr = DMPlexGetChart(plex, &pStart, &pEnd);CHKERRQ(ierr); 556 ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr); 557 558 if (patch->user_patches) { 559 /* compute patch->nuserIS, patch->userIS here */ 560 ierr = patch->userpatchconstructionop(pc, &patch->nuserIS, &patch->userIS, &patch->iterationSet, patch->userpatchconstructctx);CHKERRQ(ierr); 561 vStart = 0; 562 vEnd = patch->nuserIS; 563 } else if (patch->codim < 0) { /* codim unset */ 564 if (patch->dim < 0) { /* dim unset */ 565 ierr = DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd);CHKERRQ(ierr); 566 } else { /* dim set */ 567 ierr = DMPlexGetDepthStratum(plex, patch->dim, &vStart, &vEnd);CHKERRQ(ierr); 568 } 569 } else { /* codim set */ 570 ierr = DMPlexGetHeightStratum(plex, patch->codim, &vStart, &vEnd);CHKERRQ(ierr); 571 } 572 573 /* These labels mark the owned points. We only create patches around points that this process owns. */ 574 /* Replace this with a check of the SF */ 575 ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr); 576 ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr); 577 578 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->cellCounts);CHKERRQ(ierr); 579 cellCounts = patch->cellCounts; 580 ierr = PetscSectionSetChart(cellCounts, vStart, vEnd);CHKERRQ(ierr); 581 /* Count cells in the patch surrounding each entity */ 582 for (v = vStart; v < vEnd; ++v) { 583 PetscHashIIter hi; 584 PetscInt chtSize; 585 586 if (!patch->user_patches) { 587 ierr = DMLabelHasPoint(ghost, v, &flg);CHKERRQ(ierr); 588 /* Not an owned entity, don't make a cell patch. */ 589 if (flg) continue; 590 } 591 592 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 593 ierr = PCPatchCompleteCellPatch(dm, ht, cht);CHKERRQ(ierr); 594 PetscHashISize(cht, chtSize); 595 /* empty patch, continue */ 596 if (chtSize == 0) continue; 597 598 /* safe because size(cht) > 0 from above */ 599 PetscHashIIterBegin(cht, hi); 600 while (!PetscHashIIterAtEnd(cht, hi)) { 601 PetscInt point; 602 603 PetscHashIIterGetKey(cht, hi, point); 604 if (point >= cStart && point < cEnd) { 605 ierr = PetscSectionAddDof(cellCounts, v, 1);CHKERRQ(ierr); 606 } 607 PetscHashIIterNext(cht, hi); 608 } 609 } 610 ierr = DMLabelDestroyIndex(ghost);CHKERRQ(ierr); 611 612 ierr = PetscSectionSetUp(cellCounts);CHKERRQ(ierr); 613 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 614 ierr = PetscMalloc1(numCells, &cellsArray);CHKERRQ(ierr); 615 616 /* Now that we know how much space we need, run through again and actually remember the cells. */ 617 for (v = vStart; v < vEnd; v++ ) { 618 PetscHashIIter hi; 619 PetscInt ndof, off; 620 621 ierr = PetscSectionGetDof(cellCounts, v, &ndof);CHKERRQ(ierr); 622 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 623 if (ndof <= 0) continue; 624 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 625 ierr = PCPatchCompleteCellPatch(dm, ht, cht);CHKERRQ(ierr); 626 ndof = 0; 627 PetscHashIIterBegin(cht, hi); 628 while (!PetscHashIIterAtEnd(cht, hi)) { 629 PetscInt point; 630 631 PetscHashIIterGetKey(cht, hi, point); 632 if (point >= cStart && point < cEnd) { 633 cellsArray[ndof + off] = point; 634 ++ndof; 635 } 636 PetscHashIIterNext(cht, hi); 637 } 638 } 639 640 ierr = ISCreateGeneral(PETSC_COMM_SELF, numCells, cellsArray, PETSC_OWN_POINTER, &patch->cells);CHKERRQ(ierr); 641 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 642 patch->npatch = pEnd - pStart; 643 PetscHashIDestroy(ht); 644 PetscHashIDestroy(cht); 645 ierr = DMDestroy(&plex);CHKERRQ(ierr); 646 PetscFunctionReturn(0); 647 } 648 649 /* 650 * PCPatchCreateCellPatchDiscretisationInfo - Build the dof maps for cell patches 651 * 652 * Input Parameters: 653 * + dm - The DMPlex object defining the mesh 654 * . cellCounts - Section with counts of cells around each vertex 655 * . cells - IS of the cell point indices of cells in each patch 656 * . cellNumbering - Section mapping plex cell points to Firedrake cell indices. 657 * . nodesPerCell - number of nodes per cell. 658 * - cellNodeMap - map from cells to node indices (nodesPerCell * numCells) 659 * 660 * Output Parameters: 661 * + dofs - IS of local dof numbers of each cell in the patch 662 * . gtolCounts - Section with counts of dofs per cell patch 663 * - gtol - IS mapping from global dofs to local dofs for each patch. 664 */ 665 static PetscErrorCode PCPatchCreateCellPatchDiscretisationInfo(PC pc) 666 { 667 PC_PATCH *patch = (PC_PATCH *) pc->data; 668 PetscSection cellCounts = patch->cellCounts; 669 PetscSection gtolCounts; 670 IS cells = patch->cells; 671 PetscSection cellNumbering = patch->cellNumbering; 672 PetscInt numCells; 673 PetscInt numDofs; 674 PetscInt numGlobalDofs; 675 PetscInt totalDofsPerCell = patch->totalDofsPerCell; 676 PetscInt vStart, vEnd, v; 677 const PetscInt *cellsArray; 678 PetscInt *newCellsArray = NULL; 679 PetscInt *dofsArray = NULL; 680 PetscInt *asmArray = NULL; 681 PetscInt *globalDofsArray = NULL; 682 PetscInt globalIndex = 0; 683 PetscInt key = 0; 684 PetscInt asmKey = 0; 685 PetscHashI ht; 686 PetscErrorCode ierr; 687 688 PetscFunctionBegin; 689 /* dofcounts section is cellcounts section * dofPerCell */ 690 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 691 numDofs = numCells * totalDofsPerCell; 692 ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr); 693 ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr); 694 ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr); 695 ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr); 696 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr); 697 gtolCounts = patch->gtolCounts; 698 ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr); 699 700 ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr); 701 PetscHashICreate(ht); 702 for (v = vStart; v < vEnd; ++v) { 703 PetscInt localIndex = 0; 704 PetscInt dof, off, i, j, k, l; 705 706 PetscHashIClear(ht); 707 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 708 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 709 if (dof <= 0) continue; 710 711 for (k = 0; k < patch->nsubspaces; ++k) { 712 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 713 PetscInt nodesPerCell = patch->nodesPerCell[k]; 714 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 715 PetscInt bs = patch->bs[k]; 716 717 for (i = off; i < off + dof; ++i) { 718 /* Walk over the cells in this patch. */ 719 const PetscInt c = cellsArray[i]; 720 PetscInt cell; 721 722 ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr); 723 if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c); 724 ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr); 725 newCellsArray[i] = cell; 726 for (j = 0; j < nodesPerCell; ++j) { 727 /* For each global dof, map it into contiguous local storage. */ 728 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset; 729 /* finally, loop over block size */ 730 for (l = 0; l < bs; ++l) { 731 PetscInt localDof; 732 733 PetscHashIMap(ht, globalDof + l, localDof); 734 if (localDof == -1) { 735 localDof = localIndex++; 736 PetscHashIAdd(ht, globalDof + l, localDof); 737 } 738 if (globalIndex >= numDofs) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 739 /* And store. */ 740 dofsArray[globalIndex++] = localDof; 741 } 742 } 743 } 744 } 745 /* How many local dofs in this patch? */ 746 PetscHashISize(ht, dof); 747 ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr); 748 } 749 if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex); 750 ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr); 751 ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr); 752 ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr); 753 754 /* Now populate the global to local map. This could be merged into the above loop if we were willing to deal with reallocs. */ 755 for (v = vStart; v < vEnd; ++v) { 756 PetscHashIIter hi; 757 PetscInt dof, off, i, j, k, l; 758 759 PetscHashIClear(ht); 760 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 761 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 762 if (dof <= 0) continue; 763 764 for (k = 0; k < patch->nsubspaces; ++k) { 765 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 766 PetscInt nodesPerCell = patch->nodesPerCell[k]; 767 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 768 PetscInt bs = patch->bs[k]; 769 770 for (i = off; i < off + dof; ++i) { 771 /* Reconstruct mapping of global-to-local on this patch. */ 772 const PetscInt c = cellsArray[i]; 773 PetscInt cell; 774 775 ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr); 776 for (j = 0; j < nodesPerCell; ++j) { 777 for (l = 0; l < bs; ++l) { 778 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset + l; 779 const PetscInt localDof = dofsArray[key]; 780 781 key += 1; 782 PetscHashIAdd(ht, globalDof, localDof); 783 } 784 } 785 } 786 if (dof > 0) { 787 /* Shove it in the output data structure. */ 788 PetscInt goff; 789 790 ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr); 791 PetscHashIIterBegin(ht, hi); 792 while (!PetscHashIIterAtEnd(ht, hi)) { 793 PetscInt globalDof, localDof; 794 795 PetscHashIIterGetKeyVal(ht, hi, globalDof, localDof); 796 if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof; 797 PetscHashIIterNext(ht, hi); 798 } 799 } 800 } 801 802 /* At this point, we have a hash table ht built that maps globalDof -> localDof. 803 We need to create the dof table laid out cellwise first, then by subspace, 804 as the assembler assembles cell-wise and we need to stuff the different 805 contributions of the different function spaces to the right places. So we loop 806 over cells, then over subspaces. */ 807 if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */ 808 for (i = off; i < off + dof; ++i) { 809 const PetscInt c = cellsArray[i]; 810 PetscInt cell; 811 812 ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr); 813 for (k = 0; k < patch->nsubspaces; ++k) { 814 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 815 PetscInt nodesPerCell = patch->nodesPerCell[k]; 816 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 817 PetscInt bs = patch->bs[k]; 818 819 for (j = 0; j < nodesPerCell; ++j) { 820 for (l = 0; l < bs; ++l) { 821 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset + l; 822 PetscInt localDof; 823 824 PetscHashIMap(ht, globalDof, localDof); 825 asmArray[asmKey++] = localDof; 826 } 827 } 828 } 829 } 830 } 831 } 832 if (1 == patch->nsubspaces) {ierr = PetscMemcpy(asmArray, dofsArray, numDofs * sizeof(PetscInt));CHKERRQ(ierr);} 833 834 PetscHashIDestroy(ht); 835 ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr); 836 ierr = PetscFree(dofsArray);CHKERRQ(ierr); 837 /* Replace cell indices with firedrake-numbered ones. */ 838 ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr); 839 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr); 840 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr); 841 PetscFunctionReturn(0); 842 } 843 844 static PetscErrorCode PCPatchCreateCellPatchBCs(PC pc) 845 { 846 PC_PATCH *patch = (PC_PATCH *) pc->data; 847 const PetscInt *bcNodes = NULL; 848 PetscSection gtolCounts = patch->gtolCounts; 849 IS gtol = patch->gtol; 850 DM dm; 851 PetscInt numBcs; 852 PetscSection bcCounts; 853 PetscHashI globalBcs, localBcs, patchDofs; 854 PetscHashI ownedpts, seenpts, owneddofs, seendofs, artificialbcs; 855 PetscHashIIter hi; 856 PetscInt *bcsArray = NULL; 857 PetscInt *multBcsArray = NULL; 858 const PetscInt *gtolArray; 859 PetscInt vStart, vEnd, v, i; 860 PetscErrorCode ierr; 861 862 PetscFunctionBegin; 863 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 864 PetscHashICreate(globalBcs); 865 ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes);CHKERRQ(ierr); 866 ierr = ISGetSize(patch->ghostBcNodes, &numBcs);CHKERRQ(ierr); 867 for (i = 0; i < numBcs; ++i) { 868 /* these are already in concatenated numbering */ 869 PetscHashIAdd(globalBcs, bcNodes[i], 0); 870 } 871 ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes);CHKERRQ(ierr); 872 PetscHashICreate(patchDofs); 873 PetscHashICreate(localBcs); 874 PetscHashICreate(ownedpts); 875 PetscHashICreate(seenpts); 876 PetscHashICreate(owneddofs); 877 PetscHashICreate(seendofs); 878 PetscHashICreate(artificialbcs); 879 880 ierr = PetscSectionGetChart(patch->cellCounts, &vStart, &vEnd);CHKERRQ(ierr); 881 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->bcCounts);CHKERRQ(ierr); 882 bcCounts = patch->bcCounts; 883 ierr = PetscSectionSetChart(bcCounts, vStart, vEnd);CHKERRQ(ierr); 884 ierr = PetscMalloc1(vEnd - vStart, &patch->bcs);CHKERRQ(ierr); 885 if (patch->multiplicative) {ierr = PetscMalloc1(vEnd - vStart, &patch->multBcs);CHKERRQ(ierr);} 886 887 ierr = ISGetIndices(gtol, >olArray);CHKERRQ(ierr); 888 for (v = vStart; v < vEnd; ++v) { 889 PetscInt bcIndex = 0; 890 PetscInt multBcIndex = 0; 891 PetscInt numBcs, dof, off; 892 893 PetscHashIClear(patchDofs); 894 PetscHashIClear(localBcs); 895 ierr = PetscSectionGetDof(gtolCounts, v, &dof);CHKERRQ(ierr); 896 ierr = PetscSectionGetOffset(gtolCounts, v, &off);CHKERRQ(ierr); 897 898 if (dof <= 0) { 899 patch->bcs[v-vStart] = NULL; 900 if (patch->multiplicative) patch->multBcs[v-vStart] = NULL; 901 continue; 902 } 903 904 for (i = off; i < off + dof; ++i) { 905 const PetscInt globalDof = gtolArray[i]; 906 const PetscInt localDof = i-off; 907 PetscBool flg; 908 909 PetscHashIAdd(patchDofs, globalDof, localDof); 910 PetscHashIHasKey(globalBcs, globalDof, flg); 911 if (flg) {PetscHashIAdd(localBcs, localDof, 0);} 912 } 913 914 /* If we're doing multiplicative, make the BC data structures now corresponding solely to actual globally imposed Dirichlet BCs */ 915 if (patch->multiplicative) { 916 PetscHashISize(localBcs, numBcs); 917 ierr = PetscMalloc1(numBcs, &multBcsArray);CHKERRQ(ierr); 918 PetscHashIGetKeys(localBcs, &multBcIndex, multBcsArray); 919 ierr = PetscSortInt(numBcs, multBcsArray);CHKERRQ(ierr); 920 ierr = ISCreateGeneral(PETSC_COMM_SELF, numBcs, multBcsArray, PETSC_OWN_POINTER, &(patch->multBcs[v-vStart]));CHKERRQ(ierr); 921 } 922 923 /* Now figure out the artificial BCs: the set difference of {dofs on entities I see on the patch}\{dofs I am responsible for updating} */ 924 ierr = patch->patchconstructop((void*) patch, dm, v, ownedpts);CHKERRQ(ierr); 925 ierr = PCPatchCompleteCellPatch(dm, ownedpts, seenpts);CHKERRQ(ierr); 926 ierr = PCPatchGetPointDofs(patch, ownedpts, owneddofs, v, patch->exclude_subspace);CHKERRQ(ierr); 927 ierr = PCPatchGetPointDofs(patch, seenpts, seendofs, v, -1);CHKERRQ(ierr); 928 ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs);CHKERRQ(ierr); 929 930 if (patch->print_patches) { 931 PetscHashI globalbcdofs; 932 MPI_Comm comm; 933 934 PetscHashICreate(globalbcdofs); 935 936 ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr); 937 ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v);CHKERRQ(ierr); 938 PetscHashIIterBegin(owneddofs, hi); 939 while (!PetscHashIIterAtEnd(owneddofs, hi)) { 940 PetscInt globalDof; 941 942 PetscHashIIterGetKey(owneddofs, hi, globalDof); 943 PetscHashIIterNext(owneddofs, hi); 944 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr); 945 } 946 ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr); 947 ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v);CHKERRQ(ierr); 948 PetscHashIIterBegin(seendofs, hi); 949 while (!PetscHashIIterAtEnd(seendofs, hi)) { 950 PetscInt globalDof; 951 PetscBool flg; 952 953 PetscHashIIterGetKey(seendofs, hi, globalDof); 954 PetscHashIIterNext(seendofs, hi); 955 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr); 956 PetscHashIHasKey(globalBcs, globalDof, flg); 957 if (flg) {PetscHashIAdd(globalbcdofs, globalDof, 0);} 958 } 959 ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr); 960 ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v);CHKERRQ(ierr); 961 PetscHashISize(globalbcdofs, numBcs); 962 if (numBcs > 0) { 963 PetscHashIIterBegin(globalbcdofs, hi); 964 while (!PetscHashIIterAtEnd(globalbcdofs, hi)) { 965 PetscInt globalDof; 966 PetscHashIIterGetKey(globalbcdofs, hi, globalDof); 967 PetscHashIIterNext(globalbcdofs, hi); 968 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr); 969 } 970 } 971 ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr); 972 ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr); 973 PetscHashISize(artificialbcs, numBcs); 974 if (numBcs > 0) { 975 PetscHashIIterBegin(artificialbcs, hi); 976 while (!PetscHashIIterAtEnd(artificialbcs, hi)) { 977 PetscInt globalDof; 978 PetscHashIIterGetKey(artificialbcs, hi, globalDof); 979 PetscHashIIterNext(artificialbcs, hi); 980 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr); 981 } 982 } 983 ierr = PetscSynchronizedPrintf(comm, "\n\n");CHKERRQ(ierr); 984 ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr); 985 PetscHashIDestroy(globalbcdofs); 986 } 987 988 PetscHashISize(artificialbcs, numBcs); 989 if (numBcs > 0) { 990 PetscHashIIterBegin(artificialbcs, hi); 991 while (!PetscHashIIterAtEnd(artificialbcs, hi)) { 992 PetscInt globalDof, localDof; 993 PetscHashIIterGetKey(artificialbcs, hi, globalDof); 994 PetscHashIIterNext(artificialbcs, hi); 995 PetscHashIMap(patchDofs, globalDof, localDof); 996 if (localDof == -1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Patch %d Didn't find dof %d in patch\n", v-vStart, globalDof); 997 PetscHashIAdd(localBcs, localDof, 0); 998 } 999 } 1000 1001 /* OK, now we have a hash table with all the bcs indicated by the artificial and global bcs */ 1002 PetscHashISize(localBcs, numBcs); 1003 ierr = PetscSectionSetDof(bcCounts, v, numBcs);CHKERRQ(ierr); 1004 ierr = PetscMalloc1(numBcs, &bcsArray);CHKERRQ(ierr); 1005 PetscHashIGetKeys(localBcs, &bcIndex, bcsArray); 1006 ierr = PetscSortInt(numBcs, bcsArray);CHKERRQ(ierr); 1007 ierr = ISCreateGeneral(PETSC_COMM_SELF, numBcs, bcsArray, PETSC_OWN_POINTER, &(patch->bcs[v - vStart]));CHKERRQ(ierr); 1008 } 1009 ierr = ISRestoreIndices(gtol, >olArray);CHKERRQ(ierr); 1010 PetscHashIDestroy(artificialbcs); 1011 PetscHashIDestroy(seendofs); 1012 PetscHashIDestroy(owneddofs); 1013 PetscHashIDestroy(seenpts); 1014 PetscHashIDestroy(ownedpts); 1015 PetscHashIDestroy(localBcs); 1016 PetscHashIDestroy(patchDofs); 1017 PetscHashIDestroy(globalBcs); 1018 ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr); 1019 ierr = PetscSectionSetUp(bcCounts);CHKERRQ(ierr); 1020 PetscFunctionReturn(0); 1021 } 1022 1023 static PetscErrorCode PCPatchZeroFillMatrix_Private(Mat mat, const PetscInt ncell, const PetscInt ndof, const PetscInt *dof) 1024 { 1025 const PetscScalar *values = NULL; 1026 PetscInt rows, c, i; 1027 PetscErrorCode ierr; 1028 1029 PetscFunctionBegin; 1030 ierr = PetscCalloc1(ndof*ndof, &values);CHKERRQ(ierr); 1031 for (c = 0; c < ncell; ++c) { 1032 const PetscInt *idx = &dof[ndof*c]; 1033 ierr = MatSetValues(mat, ndof, idx, ndof, idx, values, INSERT_VALUES);CHKERRQ(ierr); 1034 } 1035 ierr = MatGetLocalSize(mat, &rows, NULL);CHKERRQ(ierr); 1036 for (i = 0; i < rows; ++i) { 1037 ierr = MatSetValues(mat, 1, &i, 1, &i, values, INSERT_VALUES);CHKERRQ(ierr); 1038 } 1039 ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1040 ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1041 ierr = PetscFree(values);CHKERRQ(ierr); 1042 PetscFunctionReturn(0); 1043 } 1044 1045 static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat) 1046 { 1047 PC_PATCH *patch = (PC_PATCH *) pc->data; 1048 Vec x, y; 1049 PetscBool flg; 1050 PetscInt csize, rsize; 1051 const char *prefix = NULL; 1052 PetscErrorCode ierr; 1053 1054 PetscFunctionBegin; 1055 x = patch->patchX[point]; 1056 y = patch->patchY[point]; 1057 ierr = VecGetSize(x, &csize);CHKERRQ(ierr); 1058 ierr = VecGetSize(y, &rsize);CHKERRQ(ierr); 1059 ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr); 1060 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 1061 ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr); 1062 ierr = MatAppendOptionsPrefix(*mat, "sub_");CHKERRQ(ierr); 1063 if (patch->sub_mat_type) {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);} 1064 ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr); 1065 ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr); 1066 if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);} 1067 /* Sparse patch matrices */ 1068 if (!flg) { 1069 PetscBT bt; 1070 PetscInt *dnnz = NULL; 1071 const PetscInt *dofsArray = NULL; 1072 PetscInt pStart, pEnd, ncell, offset, c, i, j; 1073 1074 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1075 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1076 point += pStart; 1077 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1078 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1079 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1080 ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr); 1081 ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1082 /* XXX: This uses N^2 bits to store the sparsity pattern on a 1083 * patch. This is probably OK if the patches are not too big, 1084 * but could use quite a bit of memory for planes in 3D. 1085 * Should we switch based on the value of rsize to a 1086 * hash-table (slower, but more memory efficient) approach? */ 1087 ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr); 1088 for (c = 0; c < ncell; ++c) { 1089 const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell; 1090 for (i = 0; i < patch->totalDofsPerCell; ++i) { 1091 const PetscInt row = idx[i]; 1092 for (j = 0; j < patch->totalDofsPerCell; ++j) { 1093 const PetscInt col = idx[j]; 1094 const PetscInt key = row*rsize + col; 1095 if (!PetscBTLookupSet(bt, key)) ++dnnz[row]; 1096 } 1097 } 1098 } 1099 ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 1100 ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr); 1101 ierr = PetscFree(dnnz);CHKERRQ(ierr); 1102 ierr = PCPatchZeroFillMatrix_Private(*mat, ncell, patch->totalDofsPerCell, &dofsArray[offset*patch->totalDofsPerCell]);CHKERRQ(ierr); 1103 ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1104 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1105 } 1106 ierr = MatSetUp(*mat);CHKERRQ(ierr); 1107 PetscFunctionReturn(0); 1108 } 1109 1110 static PetscErrorCode PCPatchComputeOperator_Private(PC pc, Mat mat, Mat multMat, PetscInt point) 1111 { 1112 PC_PATCH *patch = (PC_PATCH *) pc->data; 1113 const PetscInt *dofsArray; 1114 const PetscInt *cellsArray; 1115 PetscInt ncell, offset, pStart, pEnd; 1116 PetscErrorCode ierr; 1117 1118 PetscFunctionBegin; 1119 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1120 if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n"); 1121 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1122 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1123 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1124 1125 point += pStart; 1126 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1127 1128 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1129 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1130 if (ncell <= 0) { 1131 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1132 PetscFunctionReturn(0); 1133 } 1134 PetscStackPush("PCPatch user callback"); 1135 ierr = patch->usercomputeop(pc, mat, ncell, cellsArray + offset, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputectx);CHKERRQ(ierr); 1136 PetscStackPop; 1137 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1138 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1139 /* Apply boundary conditions. Could also do this through the local_to_patch guy. */ 1140 if (patch->multiplicative) { 1141 ierr = MatCopy(mat, multMat, SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1142 ierr = MatZeroRowsColumnsIS(multMat, patch->multBcs[point-pStart], (PetscScalar) 1.0, NULL, NULL);CHKERRQ(ierr); 1143 } 1144 ierr = MatZeroRowsColumnsIS(mat, patch->bcs[point-pStart], (PetscScalar) 1.0, NULL, NULL);CHKERRQ(ierr); 1145 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1146 PetscFunctionReturn(0); 1147 } 1148 1149 static PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat) 1150 { 1151 PC_PATCH *patch = (PC_PATCH *) pc->data; 1152 const PetscScalar *xArray = NULL; 1153 PetscScalar *yArray = NULL; 1154 const PetscInt *gtolArray = NULL; 1155 PetscInt dof, offset, lidx; 1156 PetscErrorCode ierr; 1157 1158 PetscFunctionBeginHot; 1159 ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 1160 ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr); 1161 ierr = VecGetArray(y, &yArray);CHKERRQ(ierr); 1162 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 1163 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 1164 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 1165 if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n"); 1166 if (mode == ADD_VALUES && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n"); 1167 for (lidx = 0; lidx < dof; ++lidx) { 1168 const PetscInt gidx = gtolArray[offset+lidx]; 1169 1170 if (mode == INSERT_VALUES) yArray[lidx] = xArray[gidx]; /* Forward */ 1171 else yArray[gidx] += xArray[lidx]; /* Reverse */ 1172 } 1173 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 1174 ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr); 1175 ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr); 1176 ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 1177 PetscFunctionReturn(0); 1178 } 1179 1180 static PetscErrorCode PCSetUp_PATCH(PC pc) 1181 { 1182 PC_PATCH *patch = (PC_PATCH *) pc->data; 1183 PetscScalar *patchX = NULL; 1184 const PetscInt *bcNodes = NULL; 1185 PetscInt numBcs, i, j; 1186 const char *prefix; 1187 PetscErrorCode ierr; 1188 1189 PetscFunctionBegin; 1190 if (!pc->setupcalled) { 1191 PetscInt pStart, pEnd, p; 1192 PetscInt localSize; 1193 1194 ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 1195 1196 localSize = patch->subspaceOffsets[patch->nsubspaces]; 1197 ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localX);CHKERRQ(ierr); 1198 ierr = VecSetUp(patch->localX);CHKERRQ(ierr); 1199 ierr = VecDuplicate(patch->localX, &patch->localY);CHKERRQ(ierr); 1200 ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr); 1201 ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr); 1202 ierr = PCPatchCreateCellPatchBCs(pc);CHKERRQ(ierr); 1203 1204 /* OK, now build the work vectors */ 1205 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr); 1206 ierr = PetscMalloc1(patch->npatch, &patch->patchX);CHKERRQ(ierr); 1207 ierr = PetscMalloc1(patch->npatch, &patch->patchY);CHKERRQ(ierr); 1208 if (patch->partition_of_unity && patch->multiplicative) {ierr = PetscMalloc1(patch->npatch, &patch->patch_dof_weights);CHKERRQ(ierr);} 1209 for (p = pStart; p < pEnd; ++p) { 1210 PetscInt dof; 1211 1212 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 1213 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchX[p-pStart]);CHKERRQ(ierr); 1214 ierr = VecSetUp(patch->patchX[p-pStart]);CHKERRQ(ierr); 1215 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchY[p-pStart]);CHKERRQ(ierr); 1216 ierr = VecSetUp(patch->patchY[p-pStart]);CHKERRQ(ierr); 1217 if (patch->partition_of_unity && patch->multiplicative) { 1218 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patch_dof_weights[p-pStart]);CHKERRQ(ierr); 1219 ierr = VecSetUp(patch->patch_dof_weights[p-pStart]);CHKERRQ(ierr); 1220 } 1221 } 1222 ierr = PetscMalloc1(patch->npatch, &patch->ksp);CHKERRQ(ierr); 1223 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 1224 for (i = 0; i < patch->npatch; ++i) { 1225 ierr = KSPCreate(PETSC_COMM_SELF, &patch->ksp[i]);CHKERRQ(ierr); 1226 ierr = KSPSetOptionsPrefix(patch->ksp[i], prefix);CHKERRQ(ierr); 1227 ierr = KSPAppendOptionsPrefix(patch->ksp[i], "sub_");CHKERRQ(ierr); 1228 } 1229 if (patch->save_operators) { 1230 ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr); 1231 if (patch->multiplicative) {ierr = PetscMalloc1(patch->npatch, &patch->multMat);CHKERRQ(ierr);} 1232 for (i = 0; i < patch->npatch; ++i) { 1233 ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i]);CHKERRQ(ierr); 1234 if (patch->multiplicative) {ierr = MatDuplicate(patch->mat[i], MAT_SHARE_NONZERO_PATTERN, &patch->multMat[i]);CHKERRQ(ierr);} 1235 } 1236 } 1237 ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 1238 1239 /* If desired, calculate weights for dof multiplicity */ 1240 if (patch->partition_of_unity) { 1241 ierr = VecDuplicate(patch->localX, &patch->dof_weights);CHKERRQ(ierr); 1242 for (i = 0; i < patch->npatch; ++i) { 1243 PetscInt dof; 1244 1245 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr); 1246 if (dof <= 0) continue; 1247 ierr = VecSet(patch->patchX[i], 1.0);CHKERRQ(ierr); 1248 /* TODO: Do we need different scatters for X and Y? */ 1249 ierr = VecGetArray(patch->patchX[i], &patchX);CHKERRQ(ierr); 1250 /* Apply bcs to patchX (zero entries) */ 1251 ierr = ISGetLocalSize(patch->bcs[i], &numBcs);CHKERRQ(ierr); 1252 ierr = ISGetIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr); 1253 for (j = 0; j < numBcs; ++j) patchX[bcNodes[j]] = 0; 1254 ierr = ISRestoreIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr); 1255 ierr = VecRestoreArray(patch->patchX[i], &patchX);CHKERRQ(ierr); 1256 1257 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchX[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 1258 } 1259 ierr = VecReciprocal(patch->dof_weights);CHKERRQ(ierr); 1260 if (patch->partition_of_unity && patch->multiplicative) { 1261 for (i = 0; i < patch->npatch; ++i) { 1262 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->dof_weights, patch->patch_dof_weights[i], INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 1263 } 1264 } 1265 } 1266 } 1267 if (patch->save_operators) { 1268 for (i = 0; i < patch->npatch; ++i) { 1269 ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr); 1270 if (patch->multiplicative) {ierr = PCPatchComputeOperator_Private(pc, patch->mat[i], patch->multMat[i], i);CHKERRQ(ierr);} 1271 else {ierr = PCPatchComputeOperator_Private(pc, patch->mat[i], NULL, i);CHKERRQ(ierr);} 1272 ierr = KSPSetOperators(patch->ksp[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr); 1273 } 1274 } 1275 if (!pc->setupcalled) { 1276 /* TODO: Only call if SetFromOptions was called on this PC */ 1277 for (i = 0; i < patch->npatch; ++i) {ierr = KSPSetFromOptions(patch->ksp[i]);CHKERRQ(ierr);} 1278 } 1279 PetscFunctionReturn(0); 1280 } 1281 1282 static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y) 1283 { 1284 PC_PATCH *patch = (PC_PATCH *) pc->data; 1285 const PetscScalar *globalX = NULL; 1286 PetscScalar *localX = NULL; 1287 PetscScalar *globalY = NULL; 1288 PetscScalar *patchX = NULL; 1289 const PetscInt *bcNodes = NULL; 1290 PetscInt nsweep = patch->symmetrise_sweep ? 2 : 1; 1291 PetscInt start[2] = {0, 0}; 1292 PetscInt end[2] = {-1, -1}; 1293 const PetscInt inc[2] = {1, -1}; 1294 const PetscScalar *localY; 1295 const PetscInt *iterationSet; 1296 PetscInt pStart, numBcs, n, sweep, bc, j; 1297 PetscErrorCode ierr; 1298 1299 PetscFunctionBegin; 1300 ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 1301 ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr); 1302 end[0] = patch->npatch; 1303 start[1] = patch->npatch-1; 1304 if (patch->user_patches) { 1305 ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr); 1306 start[1] = end[0] - 1; 1307 ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr); 1308 } 1309 /* Scatter from global space into overlapped local spaces */ 1310 ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr); 1311 ierr = VecGetArray(patch->localX, &localX);CHKERRQ(ierr); 1312 ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr); 1313 ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr); 1314 ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr); 1315 ierr = VecRestoreArray(patch->localX, &localX);CHKERRQ(ierr); 1316 1317 ierr = VecSet(patch->localY, 0.0);CHKERRQ(ierr); 1318 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr); 1319 for (sweep = 0; sweep < nsweep; sweep++) { 1320 for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) { 1321 Mat multMat = NULL; 1322 PetscInt i = patch->user_patches ? iterationSet[j] : j; 1323 PetscInt start, len; 1324 1325 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr); 1326 ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr); 1327 /* TODO: Squash out these guys in the setup as well. */ 1328 if (len <= 0) continue; 1329 /* TODO: Do we need different scatters for X and Y? */ 1330 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localX, patch->patchX[i], INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 1331 /* Apply bcs to patchX (zero entries) */ 1332 ierr = VecGetArray(patch->patchX[i], &patchX);CHKERRQ(ierr); 1333 ierr = ISGetLocalSize(patch->bcs[i], &numBcs);CHKERRQ(ierr); 1334 ierr = ISGetIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr); 1335 for (bc = 0; bc < numBcs; ++bc) patchX[bcNodes[bc]] = 0; 1336 ierr = ISRestoreIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr); 1337 ierr = VecRestoreArray(patch->patchX[i], &patchX);CHKERRQ(ierr); 1338 if (!patch->save_operators) { 1339 Mat mat; 1340 1341 ierr = PCPatchCreateMatrix_Private(pc, i, &mat);CHKERRQ(ierr); 1342 if (patch->multiplicative) {ierr = MatDuplicate(mat, MAT_SHARE_NONZERO_PATTERN, &multMat);CHKERRQ(ierr);} 1343 /* Populate operator here. */ 1344 ierr = PCPatchComputeOperator_Private(pc, mat, multMat, i);CHKERRQ(ierr); 1345 ierr = KSPSetOperators(patch->ksp[i], mat, mat); 1346 /* Drop reference so the KSPSetOperators below will blow it away. */ 1347 ierr = MatDestroy(&mat);CHKERRQ(ierr); 1348 } else if (patch->multiplicative) { 1349 multMat = patch->multMat[i]; 1350 } 1351 1352 ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 1353 ierr = KSPSolve(patch->ksp[i], patch->patchX[i], patch->patchY[i]);CHKERRQ(ierr); 1354 ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 1355 1356 if (!patch->save_operators) { 1357 PC pc; 1358 ierr = KSPSetOperators(patch->ksp[i], NULL, NULL);CHKERRQ(ierr); 1359 ierr = KSPGetPC(patch->ksp[i], &pc);CHKERRQ(ierr); 1360 /* Destroy PC context too, otherwise the factored matrix hangs around. */ 1361 ierr = PCReset(pc);CHKERRQ(ierr); 1362 } 1363 1364 if (patch->partition_of_unity && patch->multiplicative) { 1365 ierr = VecPointwiseMult(patch->patchY[i], patch->patchY[i], patch->patch_dof_weights[i]);CHKERRQ(ierr); 1366 } 1367 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchY[i], patch->localY, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 1368 1369 /* Get the matrix on the patch but with only global bcs applied. 1370 * This matrix is then multiplied with the result from the previous solve 1371 * to obtain by how much the residual changes. */ 1372 if (patch->multiplicative) { 1373 ierr = MatMult(multMat, patch->patchY[i], patch->patchX[i]);CHKERRQ(ierr); 1374 ierr = VecScale(patch->patchX[i], -1.0);CHKERRQ(ierr); 1375 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchX[i], patch->localX, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 1376 if (!patch->save_operators) {ierr = MatDestroy(&multMat);CHKERRQ(ierr);} 1377 } 1378 } 1379 } 1380 if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);} 1381 /* XXX: should we do this on the global vector? */ 1382 if (patch->partition_of_unity && !patch->multiplicative) { 1383 ierr = VecPointwiseMult(patch->localY, patch->localY, patch->dof_weights);CHKERRQ(ierr); 1384 } 1385 /* Now patch->localY contains the solution of the patch solves, so we need to combine them all. */ 1386 ierr = VecSet(y, 0.0);CHKERRQ(ierr); 1387 ierr = VecGetArray(y, &globalY);CHKERRQ(ierr); 1388 ierr = VecGetArrayRead(patch->localY, &localY);CHKERRQ(ierr); 1389 ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr); 1390 ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr); 1391 ierr = VecRestoreArrayRead(patch->localY, &localY);CHKERRQ(ierr); 1392 1393 /* Now we need to send the global BC values through */ 1394 ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr); 1395 ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr); 1396 ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 1397 ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr); 1398 for (bc = 0; bc < numBcs; ++bc) { 1399 const PetscInt idx = bcNodes[bc]; 1400 if (idx < n) globalY[idx] = globalX[idx]; 1401 } 1402 1403 ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 1404 ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr); 1405 ierr = VecRestoreArray(y, &globalY);CHKERRQ(ierr); 1406 1407 ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr); 1408 ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 1409 PetscFunctionReturn(0); 1410 } 1411 1412 static PetscErrorCode PCReset_PATCH(PC pc) 1413 { 1414 PC_PATCH *patch = (PC_PATCH *) pc->data; 1415 PetscInt i; 1416 PetscErrorCode ierr; 1417 1418 PetscFunctionBegin; 1419 /* TODO: Get rid of all these ifs */ 1420 ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr); 1421 ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr); 1422 ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr); 1423 ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr); 1424 ierr = PetscSectionDestroy(&patch->bcCounts);CHKERRQ(ierr); 1425 ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr); 1426 ierr = ISDestroy(&patch->cells);CHKERRQ(ierr); 1427 ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr); 1428 ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr); 1429 ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr); 1430 1431 if (patch->dofSection) { 1432 for (i = 0; i < patch->nsubspaces; i++) { 1433 ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr); 1434 } 1435 } 1436 ierr = PetscFree(patch->dofSection);CHKERRQ(ierr); 1437 ierr = PetscFree(patch->bs);CHKERRQ(ierr); 1438 ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr); 1439 ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr); 1440 ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr); 1441 1442 if (patch->bcs) { 1443 for (i = 0; i < patch->npatch; ++i) { 1444 ierr = ISDestroy(&patch->bcs[i]);CHKERRQ(ierr); 1445 } 1446 ierr = PetscFree(patch->bcs);CHKERRQ(ierr); 1447 } 1448 if (patch->multBcs) { 1449 for (i = 0; i < patch->npatch; ++i) { 1450 ierr = ISDestroy(&patch->multBcs[i]);CHKERRQ(ierr); 1451 } 1452 ierr = PetscFree(patch->multBcs);CHKERRQ(ierr); 1453 } 1454 if (patch->ksp) { 1455 for (i = 0; i < patch->npatch; ++i) { 1456 ierr = KSPReset(patch->ksp[i]);CHKERRQ(ierr); 1457 } 1458 } 1459 1460 ierr = VecDestroy(&patch->localX);CHKERRQ(ierr); 1461 ierr = VecDestroy(&patch->localY);CHKERRQ(ierr); 1462 if (patch->patchX) { 1463 for (i = 0; i < patch->npatch; ++i) { 1464 ierr = VecDestroy(&patch->patchX[i]);CHKERRQ(ierr); 1465 } 1466 ierr = PetscFree(patch->patchX);CHKERRQ(ierr); 1467 } 1468 if (patch->patchY) { 1469 for (i = 0; i < patch->npatch; ++i) { 1470 ierr = VecDestroy(&patch->patchY[i]);CHKERRQ(ierr); 1471 } 1472 ierr = PetscFree(patch->patchY);CHKERRQ(ierr); 1473 } 1474 if (patch->partition_of_unity) { 1475 ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr); 1476 } 1477 if (patch->patch_dof_weights) { 1478 for (i = 0; i < patch->npatch; ++i) { 1479 ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr); 1480 } 1481 ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr); 1482 } 1483 if (patch->mat) { 1484 for (i = 0; i < patch->npatch; ++i) { 1485 ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr); 1486 if (patch->multiplicative) { 1487 ierr = MatDestroy(&patch->multMat[i]);CHKERRQ(ierr); 1488 } 1489 } 1490 ierr = PetscFree(patch->mat);CHKERRQ(ierr); 1491 if (patch->multiplicative) { 1492 ierr = PetscFree(patch->multMat);CHKERRQ(ierr); 1493 } 1494 } 1495 ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr); 1496 1497 patch->bs = 0; 1498 patch->cellNodeMap = NULL; 1499 1500 if (patch->user_patches) { 1501 for (i = 0; i < patch->nuserIS; ++i) { 1502 ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr); 1503 } 1504 PetscFree(patch->userIS); 1505 patch->nuserIS = 0; 1506 } 1507 ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr); 1508 PetscFunctionReturn(0); 1509 } 1510 1511 static PetscErrorCode PCDestroy_PATCH(PC pc) 1512 { 1513 PC_PATCH *patch = (PC_PATCH *) pc->data; 1514 PetscInt i; 1515 PetscErrorCode ierr; 1516 1517 PetscFunctionBegin; 1518 ierr = PCReset_PATCH(pc);CHKERRQ(ierr); 1519 if (patch->ksp) { 1520 for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy(&patch->ksp[i]);CHKERRQ(ierr);} 1521 ierr = PetscFree(patch->ksp);CHKERRQ(ierr); 1522 } 1523 ierr = PetscFree(pc->data);CHKERRQ(ierr); 1524 PetscFunctionReturn(0); 1525 } 1526 1527 static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc) 1528 { 1529 PC_PATCH *patch = (PC_PATCH *) pc->data; 1530 PCPatchConstructType patchConstructionType = PC_PATCH_STAR; 1531 char sub_mat_type[256]; 1532 PetscBool flg, dimflg, codimflg; 1533 PetscErrorCode ierr; 1534 1535 PetscFunctionBegin; 1536 ierr = PetscOptionsHead(PetscOptionsObject, "Vertex-patch Additive Schwarz options");CHKERRQ(ierr); 1537 ierr = PetscOptionsBool("-pc_patch_save_operators", "Store all patch operators for lifetime of PC?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr); 1538 ierr = PetscOptionsBool("-pc_patch_partition_of_unity", "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr); 1539 ierr = PetscOptionsBool("-pc_patch_multiplicative", "Gauss-Seidel instead of Jacobi?", "PCPatchSetMultiplicative", patch->multiplicative, &patch->multiplicative, &flg);CHKERRQ(ierr); 1540 ierr = PetscOptionsInt("-pc_patch_construct_dim", "What dimension of mesh point to construct patches by? (0 = vertices)", "PCSetFromOptions_PATCH", patch->dim, &patch->dim, &dimflg);CHKERRQ(ierr); 1541 ierr = PetscOptionsInt("-pc_patch_construct_codim", "What co-dimension of mesh point to construct patches by? (0 = cells)", "PCSetFromOptions_PATCH", patch->codim, &patch->codim, &codimflg);CHKERRQ(ierr); 1542 if (dimflg && codimflg) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension");CHKERRQ(ierr); 1543 ierr = PetscOptionsEnum("-pc_patch_construct_type", "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr); 1544 ierr = PetscOptionsInt("-pc_patch_vanka_dim", "Topological dimension of entities for Vanka to ignore", "PCSetFromOptions_PATCH", patch->vankadim, &patch->vankadim, &flg); 1545 if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);} 1546 ierr = PetscOptionsFList("-pc_patch_sub_mat_type", "Matrix type for patch solves", "PCPatchSetSubMatType", MatList, NULL, sub_mat_type, 256, &flg);CHKERRQ(ierr); 1547 if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);} 1548 ierr = PetscOptionsBool("-pc_patch_print_patches", "Print out information during patch construction?", "PCSetFromOptions_PATCH", patch->print_patches, &patch->print_patches, &flg);CHKERRQ(ierr); 1549 ierr = PetscOptionsBool("-pc_patch_symmetrise_sweep", "Go start->end, end->start?", "PCSetFromOptions_PATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr); 1550 ierr = PetscOptionsInt("-pc_patch_exclude_subspace", "What subspace (if any) to exclude in construction?", "PCSetFromOptions_PATCH", patch->exclude_subspace, &patch->exclude_subspace, &flg); 1551 ierr = PetscOptionsTail();CHKERRQ(ierr); 1552 PetscFunctionReturn(0); 1553 } 1554 1555 static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc) 1556 { 1557 PC_PATCH *patch = (PC_PATCH*) pc->data; 1558 KSPConvergedReason reason; 1559 PetscInt i; 1560 PetscErrorCode ierr; 1561 1562 PetscFunctionBegin; 1563 for (i = 0; i < patch->npatch; ++i) { 1564 ierr = KSPSetUp(patch->ksp[i]);CHKERRQ(ierr); 1565 ierr = KSPGetConvergedReason(patch->ksp[i], &reason);CHKERRQ(ierr); 1566 if (reason == KSP_DIVERGED_PCSETUP_FAILED) pc->failedreason = PC_SUBPC_ERROR; 1567 } 1568 PetscFunctionReturn(0); 1569 } 1570 1571 static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer) 1572 { 1573 PC_PATCH *patch = (PC_PATCH *) pc->data; 1574 PetscViewer sviewer; 1575 PetscBool isascii; 1576 PetscMPIInt rank; 1577 PetscErrorCode ierr; 1578 1579 PetscFunctionBegin; 1580 /* TODO Redo tabbing with set tbas in new style */ 1581 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 1582 if (!isascii) PetscFunctionReturn(0); 1583 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr); 1584 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1585 ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr); 1586 if (patch->multiplicative) {ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: multiplicative\n");CHKERRQ(ierr);} 1587 else {ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr);} 1588 if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);} 1589 else {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);} 1590 if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);} 1591 else {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);} 1592 if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);} 1593 else {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);} 1594 if (patch->patchconstructop == PCPatchConstruct_Star) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);} 1595 else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);} 1596 else if (patch->patchconstructop == PCPatchConstruct_User) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);} 1597 else {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);} 1598 ierr = PetscViewerASCIIPrintf(viewer, "KSP on patches (all same):\n");CHKERRQ(ierr); 1599 if (patch->ksp) { 1600 ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 1601 if (!rank) { 1602 ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr); 1603 ierr = KSPView(patch->ksp[0], sviewer);CHKERRQ(ierr); 1604 ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr); 1605 } 1606 ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 1607 } else { 1608 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1609 ierr = PetscViewerASCIIPrintf(viewer, "KSP not yet set.\n");CHKERRQ(ierr); 1610 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1611 } 1612 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1613 PetscFunctionReturn(0); 1614 } 1615 1616 /*MC 1617 PCPATCH = "patch" - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping 1618 small block additive and multiplicative preconditioners. Block definition is based on topology from 1619 a DM and equation numbering from a PetscSection. 1620 1621 Options Database Keys: 1622 + -pc_patch_cells_view - Views the process local cell numbers for each patch 1623 . -pc_patch_points_view - Views the process local mesh point numbers for each patch 1624 . -pc_patch_g2l_view - Views the map between global dofs and patch local dofs for each patch 1625 . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary 1626 - -pc_patch_sub_mat_view - Views the matrix associated with each patch 1627 1628 Level: intermediate 1629 1630 .seealso: PCType, PCCreate(), PCSetType() 1631 M*/ 1632 PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc) 1633 { 1634 PC_PATCH *patch; 1635 PetscErrorCode ierr; 1636 1637 PetscFunctionBegin; 1638 ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr); 1639 1640 /* Set some defaults */ 1641 patch->sub_mat_type = NULL; 1642 patch->save_operators = PETSC_TRUE; 1643 patch->partition_of_unity = PETSC_FALSE; 1644 patch->multiplicative = PETSC_FALSE; 1645 patch->codim = -1; 1646 patch->dim = -1; 1647 patch->exclude_subspace = -1; 1648 patch->vankadim = -1; 1649 patch->patchconstructop = PCPatchConstruct_Star; 1650 patch->print_patches = PETSC_FALSE; 1651 patch->symmetrise_sweep = PETSC_FALSE; 1652 patch->nuserIS = 0; 1653 patch->userIS = NULL; 1654 patch->iterationSet = NULL; 1655 patch->user_patches = PETSC_FALSE; 1656 1657 pc->data = (void *) patch; 1658 pc->ops->apply = PCApply_PATCH; 1659 pc->ops->applytranspose = 0; /* PCApplyTranspose_PATCH; */ 1660 pc->ops->setup = PCSetUp_PATCH; 1661 pc->ops->reset = PCReset_PATCH; 1662 pc->ops->destroy = PCDestroy_PATCH; 1663 pc->ops->setfromoptions = PCSetFromOptions_PATCH; 1664 pc->ops->setuponblocks = PCSetUpOnBlocks_PATCH; 1665 pc->ops->view = PCView_PATCH; 1666 pc->ops->applyrichardson = 0; 1667 1668 PetscFunctionReturn(0); 1669 } 1670