1 #include <petsc/private/pcpatchimpl.h> /*I "petscpc.h" I*/ 2 #include <petsc/private/kspimpl.h> /* For ksp->setfromoptionscalled */ 3 #include <petsc/private/dmpleximpl.h> /* For DMPlexComputeJacobian_Patch_Internal() */ 4 #include <petscsf.h> 5 #include <petscbt.h> 6 #include <petscds.h> 7 8 PetscLogEvent PC_Patch_CreatePatches, PC_Patch_ComputeOp, PC_Patch_Solve, PC_Patch_Scatter, PC_Patch_Apply, PC_Patch_Prealloc; 9 10 PETSC_STATIC_INLINE PetscErrorCode ObjectView(PetscObject obj, PetscViewer viewer, PetscViewerFormat format) 11 { 12 PetscErrorCode ierr; 13 14 ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 15 ierr = PetscObjectView(obj, viewer);CHKERRQ(ierr); 16 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 17 return(0); 18 } 19 20 static PetscErrorCode PCPatchConstruct_Star(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 21 { 22 PetscInt starSize; 23 PetscInt *star = NULL, si; 24 PetscErrorCode ierr; 25 26 PetscFunctionBegin; 27 PetscHSetIClear(ht); 28 /* To start with, add the point we care about */ 29 ierr = PetscHSetIAdd(ht, point);CHKERRQ(ierr); 30 /* Loop over all the points that this point connects to */ 31 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 32 for (si = 0; si < starSize*2; si += 2) {ierr = PetscHSetIAdd(ht, star[si]);CHKERRQ(ierr);} 33 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 34 PetscFunctionReturn(0); 35 } 36 37 static PetscErrorCode PCPatchConstruct_Vanka(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 38 { 39 PC_PATCH *patch = (PC_PATCH *) vpatch; 40 PetscInt starSize; 41 PetscInt *star = NULL; 42 PetscBool shouldIgnore = PETSC_FALSE; 43 PetscInt cStart, cEnd, iStart, iEnd, si; 44 PetscErrorCode ierr; 45 46 PetscFunctionBegin; 47 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 48 /* To start with, add the point we care about */ 49 ierr = PetscHSetIAdd(ht, point);CHKERRQ(ierr); 50 /* Should we ignore any points of a certain dimension? */ 51 if (patch->vankadim >= 0) { 52 shouldIgnore = PETSC_TRUE; 53 ierr = DMPlexGetDepthStratum(dm, patch->vankadim, &iStart, &iEnd);CHKERRQ(ierr); 54 } 55 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 56 /* Loop over all the cells that this point connects to */ 57 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 58 for (si = 0; si < starSize*2; si += 2) { 59 const PetscInt cell = star[si]; 60 PetscInt closureSize; 61 PetscInt *closure = NULL, ci; 62 63 if (cell < cStart || cell >= cEnd) continue; 64 /* now loop over all entities in the closure of that cell */ 65 ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 66 for (ci = 0; ci < closureSize*2; ci += 2) { 67 const PetscInt newpoint = closure[ci]; 68 69 /* We've been told to ignore entities of this type.*/ 70 if (shouldIgnore && newpoint >= iStart && newpoint < iEnd) continue; 71 ierr = PetscHSetIAdd(ht, newpoint);CHKERRQ(ierr); 72 } 73 ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 74 } 75 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 76 PetscFunctionReturn(0); 77 } 78 79 static PetscErrorCode PCPatchConstruct_Pardecomp(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 80 { 81 PC_PATCH *patch = (PC_PATCH *) vpatch; 82 DMLabel ghost = NULL; 83 const PetscInt *leaves; 84 PetscInt nleaves, pStart, pEnd, loc; 85 PetscBool isFiredrake; 86 DM plex; 87 PetscBool flg; 88 PetscInt starSize; 89 PetscInt *star = NULL; 90 PetscInt opoint, overlapi; 91 PetscErrorCode ierr; 92 93 PetscFunctionBegin; 94 PetscHSetIClear(ht); 95 96 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 97 ierr = DMPlexGetChart(plex, &pStart, &pEnd);CHKERRQ(ierr); 98 99 ierr = DMHasLabel(dm, "pyop2_ghost", &isFiredrake);CHKERRQ(ierr); 100 if (isFiredrake) { 101 ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr); 102 ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr); 103 } else { 104 PetscSF sf; 105 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 106 ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 107 nleaves = PetscMax(nleaves, 0); 108 } 109 110 for (opoint = pStart; opoint < pEnd; ++opoint) { 111 if (ghost) {ierr = DMLabelHasPoint(ghost, opoint, &flg);CHKERRQ(ierr);} 112 else {ierr = PetscFindInt(opoint, nleaves, leaves, &loc);CHKERRQ(ierr); flg = loc >=0 ? PETSC_TRUE : PETSC_FALSE;} 113 /* Not an owned entity, don't make a cell patch. */ 114 if (flg) continue; 115 ierr = PetscHSetIAdd(ht, opoint);CHKERRQ(ierr); 116 } 117 118 /* Now build the overlap for the patch */ 119 for (overlapi = 0; overlapi < patch->pardecomp_overlap; ++overlapi) { 120 PetscInt index = 0; 121 PetscInt *htpoints = NULL; 122 PetscInt htsize; 123 PetscInt i; 124 125 ierr = PetscHSetIGetSize(ht, &htsize);CHKERRQ(ierr); 126 ierr = PetscMalloc1(htsize, &htpoints);CHKERRQ(ierr); 127 ierr = PetscHSetIGetElems(ht, &index, htpoints);CHKERRQ(ierr); 128 129 for (i = 0; i < htsize; ++i) { 130 PetscInt hpoint = htpoints[i]; 131 PetscInt si; 132 133 ierr = DMPlexGetTransitiveClosure(dm, hpoint, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 134 for (si = 0; si < starSize*2; si += 2) { 135 const PetscInt starp = star[si]; 136 PetscInt closureSize; 137 PetscInt *closure = NULL, ci; 138 139 /* now loop over all entities in the closure of starp */ 140 ierr = DMPlexGetTransitiveClosure(dm, starp, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 141 for (ci = 0; ci < closureSize*2; ci += 2) { 142 const PetscInt closstarp = closure[ci]; 143 ierr = PetscHSetIAdd(ht, closstarp);CHKERRQ(ierr); 144 } 145 ierr = DMPlexRestoreTransitiveClosure(dm, starp, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 146 } 147 ierr = DMPlexRestoreTransitiveClosure(dm, hpoint, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 148 } 149 ierr = PetscFree(htpoints);CHKERRQ(ierr); 150 } 151 152 PetscFunctionReturn(0); 153 } 154 155 /* The user's already set the patches in patch->userIS. Build the hash tables */ 156 static PetscErrorCode PCPatchConstruct_User(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 157 { 158 PC_PATCH *patch = (PC_PATCH *) vpatch; 159 IS patchis = patch->userIS[point]; 160 PetscInt n; 161 const PetscInt *patchdata; 162 PetscInt pStart, pEnd, i; 163 PetscErrorCode ierr; 164 165 PetscFunctionBegin; 166 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 167 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 168 ierr = ISGetLocalSize(patchis, &n);CHKERRQ(ierr); 169 ierr = ISGetIndices(patchis, &patchdata);CHKERRQ(ierr); 170 for (i = 0; i < n; ++i) { 171 const PetscInt ownedpoint = patchdata[i]; 172 173 if (ownedpoint < pStart || ownedpoint >= pEnd) { 174 SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D was not in [%D, %D)", ownedpoint, pStart, pEnd); 175 } 176 ierr = PetscHSetIAdd(ht, ownedpoint);CHKERRQ(ierr); 177 } 178 ierr = ISRestoreIndices(patchis, &patchdata);CHKERRQ(ierr); 179 PetscFunctionReturn(0); 180 } 181 182 static PetscErrorCode PCPatchCreateDefaultSF_Private(PC pc, PetscInt n, const PetscSF *sf, const PetscInt *bs) 183 { 184 PC_PATCH *patch = (PC_PATCH *) pc->data; 185 PetscInt i; 186 PetscErrorCode ierr; 187 188 PetscFunctionBegin; 189 if (n == 1 && bs[0] == 1) { 190 patch->defaultSF = sf[0]; 191 ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr); 192 } else { 193 PetscInt allRoots = 0, allLeaves = 0; 194 PetscInt leafOffset = 0; 195 PetscInt *ilocal = NULL; 196 PetscSFNode *iremote = NULL; 197 PetscInt *remoteOffsets = NULL; 198 PetscInt index = 0; 199 PetscHMapI rankToIndex; 200 PetscInt numRanks = 0; 201 PetscSFNode *remote = NULL; 202 PetscSF rankSF; 203 PetscInt *ranks = NULL; 204 PetscInt *offsets = NULL; 205 MPI_Datatype contig; 206 PetscHSetI ranksUniq; 207 208 /* First figure out how many dofs there are in the concatenated numbering. 209 * allRoots: number of owned global dofs; 210 * allLeaves: number of visible dofs (global + ghosted). 211 */ 212 for (i = 0; i < n; ++i) { 213 PetscInt nroots, nleaves; 214 215 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 216 allRoots += nroots * bs[i]; 217 allLeaves += nleaves * bs[i]; 218 } 219 ierr = PetscMalloc1(allLeaves, &ilocal);CHKERRQ(ierr); 220 ierr = PetscMalloc1(allLeaves, &iremote);CHKERRQ(ierr); 221 /* Now build an SF that just contains process connectivity. */ 222 ierr = PetscHSetICreate(&ranksUniq);CHKERRQ(ierr); 223 for (i = 0; i < n; ++i) { 224 const PetscMPIInt *ranks = NULL; 225 PetscInt nranks, j; 226 227 ierr = PetscSFSetUp(sf[i]);CHKERRQ(ierr); 228 ierr = PetscSFGetRanks(sf[i], &nranks, &ranks, NULL, NULL, NULL);CHKERRQ(ierr); 229 /* These are all the ranks who communicate with me. */ 230 for (j = 0; j < nranks; ++j) { 231 ierr = PetscHSetIAdd(ranksUniq, (PetscInt) ranks[j]);CHKERRQ(ierr); 232 } 233 } 234 ierr = PetscHSetIGetSize(ranksUniq, &numRanks);CHKERRQ(ierr); 235 ierr = PetscMalloc1(numRanks, &remote);CHKERRQ(ierr); 236 ierr = PetscMalloc1(numRanks, &ranks);CHKERRQ(ierr); 237 ierr = PetscHSetIGetElems(ranksUniq, &index, ranks);CHKERRQ(ierr); 238 239 ierr = PetscHMapICreate(&rankToIndex);CHKERRQ(ierr); 240 for (i = 0; i < numRanks; ++i) { 241 remote[i].rank = ranks[i]; 242 remote[i].index = 0; 243 ierr = PetscHMapISet(rankToIndex, ranks[i], i);CHKERRQ(ierr); 244 } 245 ierr = PetscFree(ranks);CHKERRQ(ierr); 246 ierr = PetscHSetIDestroy(&ranksUniq);CHKERRQ(ierr); 247 ierr = PetscSFCreate(PetscObjectComm((PetscObject) pc), &rankSF);CHKERRQ(ierr); 248 ierr = PetscSFSetGraph(rankSF, 1, numRanks, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 249 ierr = PetscSFSetUp(rankSF);CHKERRQ(ierr); 250 /* OK, use it to communicate the root offset on the remote 251 * processes for each subspace. */ 252 ierr = PetscMalloc1(n, &offsets);CHKERRQ(ierr); 253 ierr = PetscMalloc1(n*numRanks, &remoteOffsets);CHKERRQ(ierr); 254 255 offsets[0] = 0; 256 for (i = 1; i < n; ++i) { 257 PetscInt nroots; 258 259 ierr = PetscSFGetGraph(sf[i-1], &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 260 offsets[i] = offsets[i-1] + nroots*bs[i-1]; 261 } 262 /* Offsets are the offsets on the current process of the 263 * global dof numbering for the subspaces. */ 264 ierr = MPI_Type_contiguous(n, MPIU_INT, &contig);CHKERRQ(ierr); 265 ierr = MPI_Type_commit(&contig);CHKERRQ(ierr); 266 267 ierr = PetscSFBcastBegin(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 268 ierr = PetscSFBcastEnd(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 269 ierr = MPI_Type_free(&contig);CHKERRQ(ierr); 270 ierr = PetscFree(offsets);CHKERRQ(ierr); 271 ierr = PetscSFDestroy(&rankSF);CHKERRQ(ierr); 272 /* Now remoteOffsets contains the offsets on the remote 273 * processes who communicate with me. So now we can 274 * concatenate the list of SFs into a single one. */ 275 index = 0; 276 for (i = 0; i < n; ++i) { 277 const PetscSFNode *remote = NULL; 278 const PetscInt *local = NULL; 279 PetscInt nroots, nleaves, j; 280 281 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, &local, &remote);CHKERRQ(ierr); 282 for (j = 0; j < nleaves; ++j) { 283 PetscInt rank = remote[j].rank; 284 PetscInt idx, rootOffset, k; 285 286 ierr = PetscHMapIGet(rankToIndex, rank, &idx);CHKERRQ(ierr); 287 if (idx == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Didn't find rank, huh?"); 288 /* Offset on given rank for ith subspace */ 289 rootOffset = remoteOffsets[n*idx + i]; 290 for (k = 0; k < bs[i]; ++k) { 291 ilocal[index] = (local ? local[j] : j)*bs[i] + k + leafOffset; 292 iremote[index].rank = remote[j].rank; 293 iremote[index].index = remote[j].index*bs[i] + k + rootOffset; 294 ++index; 295 } 296 } 297 leafOffset += nleaves * bs[i]; 298 } 299 ierr = PetscHMapIDestroy(&rankToIndex);CHKERRQ(ierr); 300 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 301 ierr = PetscSFCreate(PetscObjectComm((PetscObject)pc), &patch->defaultSF);CHKERRQ(ierr); 302 ierr = PetscSFSetGraph(patch->defaultSF, allRoots, allLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 303 } 304 PetscFunctionReturn(0); 305 } 306 307 /* TODO: Docs */ 308 PetscErrorCode PCPatchSetIgnoreDim(PC pc, PetscInt dim) 309 { 310 PC_PATCH *patch = (PC_PATCH *) pc->data; 311 PetscFunctionBegin; 312 patch->ignoredim = dim; 313 PetscFunctionReturn(0); 314 } 315 316 /* TODO: Docs */ 317 PetscErrorCode PCPatchGetIgnoreDim(PC pc, PetscInt *dim) 318 { 319 PC_PATCH *patch = (PC_PATCH *) pc->data; 320 PetscFunctionBegin; 321 *dim = patch->ignoredim; 322 PetscFunctionReturn(0); 323 } 324 325 /* TODO: Docs */ 326 PetscErrorCode PCPatchSetSaveOperators(PC pc, PetscBool flg) 327 { 328 PC_PATCH *patch = (PC_PATCH *) pc->data; 329 PetscFunctionBegin; 330 patch->save_operators = flg; 331 PetscFunctionReturn(0); 332 } 333 334 /* TODO: Docs */ 335 PetscErrorCode PCPatchGetSaveOperators(PC pc, PetscBool *flg) 336 { 337 PC_PATCH *patch = (PC_PATCH *) pc->data; 338 PetscFunctionBegin; 339 *flg = patch->save_operators; 340 PetscFunctionReturn(0); 341 } 342 343 /* TODO: Docs */ 344 PetscErrorCode PCPatchSetPrecomputeElementTensors(PC pc, PetscBool flg) 345 { 346 PC_PATCH *patch = (PC_PATCH *) pc->data; 347 PetscFunctionBegin; 348 patch->precomputeElementTensors = flg; 349 PetscFunctionReturn(0); 350 } 351 352 /* TODO: Docs */ 353 PetscErrorCode PCPatchGetPrecomputeElementTensors(PC pc, PetscBool *flg) 354 { 355 PC_PATCH *patch = (PC_PATCH *) pc->data; 356 PetscFunctionBegin; 357 *flg = patch->precomputeElementTensors; 358 PetscFunctionReturn(0); 359 } 360 361 /* TODO: Docs */ 362 PetscErrorCode PCPatchSetPartitionOfUnity(PC pc, PetscBool flg) 363 { 364 PC_PATCH *patch = (PC_PATCH *) pc->data; 365 PetscFunctionBegin; 366 patch->partition_of_unity = flg; 367 PetscFunctionReturn(0); 368 } 369 370 /* TODO: Docs */ 371 PetscErrorCode PCPatchGetPartitionOfUnity(PC pc, PetscBool *flg) 372 { 373 PC_PATCH *patch = (PC_PATCH *) pc->data; 374 PetscFunctionBegin; 375 *flg = patch->partition_of_unity; 376 PetscFunctionReturn(0); 377 } 378 379 /* TODO: Docs */ 380 PetscErrorCode PCPatchSetLocalComposition(PC pc, PCCompositeType type) 381 { 382 PC_PATCH *patch = (PC_PATCH *) pc->data; 383 PetscFunctionBegin; 384 if (type != PC_COMPOSITE_ADDITIVE && type != PC_COMPOSITE_MULTIPLICATIVE) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Only supports additive or multiplicative as the local type"); 385 patch->local_composition_type = type; 386 PetscFunctionReturn(0); 387 } 388 389 /* TODO: Docs */ 390 PetscErrorCode PCPatchGetLocalComposition(PC pc, PCCompositeType *type) 391 { 392 PC_PATCH *patch = (PC_PATCH *) pc->data; 393 PetscFunctionBegin; 394 *type = patch->local_composition_type; 395 PetscFunctionReturn(0); 396 } 397 398 /* TODO: Docs */ 399 PetscErrorCode PCPatchSetSubMatType(PC pc, MatType sub_mat_type) 400 { 401 PC_PATCH *patch = (PC_PATCH *) pc->data; 402 PetscErrorCode ierr; 403 404 PetscFunctionBegin; 405 if (patch->sub_mat_type) {ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);} 406 ierr = PetscStrallocpy(sub_mat_type, (char **) &patch->sub_mat_type);CHKERRQ(ierr); 407 PetscFunctionReturn(0); 408 } 409 410 /* TODO: Docs */ 411 PetscErrorCode PCPatchGetSubMatType(PC pc, MatType *sub_mat_type) 412 { 413 PC_PATCH *patch = (PC_PATCH *) pc->data; 414 PetscFunctionBegin; 415 *sub_mat_type = patch->sub_mat_type; 416 PetscFunctionReturn(0); 417 } 418 419 /* TODO: Docs */ 420 PetscErrorCode PCPatchSetCellNumbering(PC pc, PetscSection cellNumbering) 421 { 422 PC_PATCH *patch = (PC_PATCH *) pc->data; 423 PetscErrorCode ierr; 424 425 PetscFunctionBegin; 426 patch->cellNumbering = cellNumbering; 427 ierr = PetscObjectReference((PetscObject) cellNumbering);CHKERRQ(ierr); 428 PetscFunctionReturn(0); 429 } 430 431 /* TODO: Docs */ 432 PetscErrorCode PCPatchGetCellNumbering(PC pc, PetscSection *cellNumbering) 433 { 434 PC_PATCH *patch = (PC_PATCH *) pc->data; 435 PetscFunctionBegin; 436 *cellNumbering = patch->cellNumbering; 437 PetscFunctionReturn(0); 438 } 439 440 /* TODO: Docs */ 441 PetscErrorCode PCPatchSetConstructType(PC pc, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx) 442 { 443 PC_PATCH *patch = (PC_PATCH *) pc->data; 444 445 PetscFunctionBegin; 446 patch->ctype = ctype; 447 switch (ctype) { 448 case PC_PATCH_STAR: 449 patch->user_patches = PETSC_FALSE; 450 patch->patchconstructop = PCPatchConstruct_Star; 451 break; 452 case PC_PATCH_VANKA: 453 patch->user_patches = PETSC_FALSE; 454 patch->patchconstructop = PCPatchConstruct_Vanka; 455 break; 456 case PC_PATCH_PARDECOMP: 457 patch->user_patches = PETSC_FALSE; 458 patch->patchconstructop = PCPatchConstruct_Pardecomp; 459 break; 460 case PC_PATCH_USER: 461 case PC_PATCH_PYTHON: 462 patch->user_patches = PETSC_TRUE; 463 patch->patchconstructop = PCPatchConstruct_User; 464 if (func) { 465 patch->userpatchconstructionop = func; 466 patch->userpatchconstructctx = ctx; 467 } 468 break; 469 default: 470 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 471 } 472 PetscFunctionReturn(0); 473 } 474 475 /* TODO: Docs */ 476 PetscErrorCode PCPatchGetConstructType(PC pc, PCPatchConstructType *ctype, PetscErrorCode (**func)(PC, PetscInt *, IS **, IS *, void *), void **ctx) 477 { 478 PC_PATCH *patch = (PC_PATCH *) pc->data; 479 480 PetscFunctionBegin; 481 *ctype = patch->ctype; 482 switch (patch->ctype) { 483 case PC_PATCH_STAR: 484 case PC_PATCH_VANKA: 485 case PC_PATCH_PARDECOMP: 486 break; 487 case PC_PATCH_USER: 488 case PC_PATCH_PYTHON: 489 *func = patch->userpatchconstructionop; 490 *ctx = patch->userpatchconstructctx; 491 break; 492 default: 493 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 494 } 495 PetscFunctionReturn(0); 496 } 497 498 /* TODO: Docs */ 499 PetscErrorCode PCPatchSetDiscretisationInfo(PC pc, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, 500 const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes) 501 { 502 PC_PATCH *patch = (PC_PATCH *) pc->data; 503 DM dm; 504 PetscSF *sfs; 505 PetscInt cStart, cEnd, i, j; 506 PetscErrorCode ierr; 507 508 PetscFunctionBegin; 509 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 510 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 511 ierr = PetscMalloc1(nsubspaces, &sfs);CHKERRQ(ierr); 512 ierr = PetscMalloc1(nsubspaces, &patch->dofSection);CHKERRQ(ierr); 513 ierr = PetscMalloc1(nsubspaces, &patch->bs);CHKERRQ(ierr); 514 ierr = PetscMalloc1(nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr); 515 ierr = PetscMalloc1(nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr); 516 ierr = PetscMalloc1(nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr); 517 518 patch->nsubspaces = nsubspaces; 519 patch->totalDofsPerCell = 0; 520 for (i = 0; i < nsubspaces; ++i) { 521 ierr = DMGetDefaultSection(dms[i], &patch->dofSection[i]);CHKERRQ(ierr); 522 ierr = PetscObjectReference((PetscObject) patch->dofSection[i]);CHKERRQ(ierr); 523 ierr = DMGetDefaultSF(dms[i], &sfs[i]);CHKERRQ(ierr); 524 patch->bs[i] = bs[i]; 525 patch->nodesPerCell[i] = nodesPerCell[i]; 526 patch->totalDofsPerCell += nodesPerCell[i]*bs[i]; 527 ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr); 528 for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j]; 529 patch->subspaceOffsets[i] = subspaceOffsets[i]; 530 } 531 ierr = PCPatchCreateDefaultSF_Private(pc, nsubspaces, sfs, patch->bs);CHKERRQ(ierr); 532 ierr = PetscFree(sfs);CHKERRQ(ierr); 533 534 patch->subspaceOffsets[nsubspaces] = subspaceOffsets[nsubspaces]; 535 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr); 536 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr); 537 PetscFunctionReturn(0); 538 } 539 540 /* TODO: Docs */ 541 PetscErrorCode PCPatchSetDiscretisationInfoCombined(PC pc, DM dm, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes) 542 { 543 PC_PATCH *patch = (PC_PATCH *) pc->data; 544 PetscInt cStart, cEnd, i, j; 545 PetscErrorCode ierr; 546 547 PetscFunctionBegin; 548 patch->combined = PETSC_TRUE; 549 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 550 ierr = DMGetNumFields(dm, &patch->nsubspaces);CHKERRQ(ierr); 551 ierr = PetscCalloc1(patch->nsubspaces, &patch->dofSection);CHKERRQ(ierr); 552 ierr = PetscMalloc1(patch->nsubspaces, &patch->bs);CHKERRQ(ierr); 553 ierr = PetscMalloc1(patch->nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr); 554 ierr = PetscMalloc1(patch->nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr); 555 ierr = PetscCalloc1(patch->nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr); 556 ierr = DMGetDefaultSection(dm, &patch->dofSection[0]);CHKERRQ(ierr); 557 ierr = PetscObjectReference((PetscObject) patch->dofSection[0]);CHKERRQ(ierr); 558 ierr = PetscSectionGetStorageSize(patch->dofSection[0], &patch->subspaceOffsets[patch->nsubspaces]);CHKERRQ(ierr); 559 patch->totalDofsPerCell = 0; 560 for (i = 0; i < patch->nsubspaces; ++i) { 561 patch->bs[i] = 1; 562 patch->nodesPerCell[i] = nodesPerCell[i]; 563 patch->totalDofsPerCell += nodesPerCell[i]; 564 ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr); 565 for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j]; 566 } 567 ierr = DMGetDefaultSF(dm, &patch->defaultSF);CHKERRQ(ierr); 568 ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr); 569 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr); 570 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr); 571 PetscFunctionReturn(0); 572 } 573 574 /*@C 575 576 PCPatchSetComputeFunction - Set the callback used to compute patch residuals 577 578 Input Parameters: 579 + pc - The PC 580 . func - The callback 581 - ctx - The user context 582 583 Calling sequence of func: 584 $ func (PC pc,PetscInt point,Vec x,Vec f,IS cellIS,PetscInt n,const PetscInt* dofsArray,const PetscInt* dofsArrayWithAll,void* ctx) 585 586 + pc - The PC 587 . point - The point 588 . x - The input solution (not used in linear problems) 589 . f - The patch residual vector 590 . cellIS - An array of the cell numbers 591 . n - The size of dofsArray 592 . dofsArray - The dofmap for the dofs to be solved for 593 . dofsArrayWithAll - The dofmap for all dofs on the patch 594 - ctx - The user context 595 596 Level: advanced 597 598 Notes: 599 The matrix entries have been set to zero before the call. 600 601 .seealso: PCPatchSetComputeOperator(), PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo() 602 @*/ 603 PetscErrorCode PCPatchSetComputeFunction(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 604 { 605 PC_PATCH *patch = (PC_PATCH *) pc->data; 606 607 PetscFunctionBegin; 608 patch->usercomputef = func; 609 patch->usercomputefctx = ctx; 610 PetscFunctionReturn(0); 611 } 612 613 /*@C 614 615 PCPatchSetComputeFunctionInteriorFacets - Set the callback used to compute facet integrals for patch residuals 616 617 Logically collective on PC 618 619 Input Parameters: 620 + pc - The PC 621 . func - The callback 622 - ctx - The user context 623 624 Calling sequence of func: 625 $ func (PC pc,PetscInt point,Vec x,Vec f,IS facetIS,PetscInt n,const PetscInt* dofsArray,const PetscInt* dofsArrayWithAll,void* ctx) 626 627 + pc - The PC 628 . point - The point 629 . x - The input solution (not used in linear problems) 630 . f - The patch residual vector 631 . facetIS - An array of the facet numbers 632 . n - The size of dofsArray 633 . dofsArray - The dofmap for the dofs to be solved for 634 . dofsArrayWithAll - The dofmap for all dofs on the patch 635 - ctx - The user context 636 637 Level: advanced 638 639 Notes: 640 The matrix entries have been set to zero before the call. 641 642 .seealso: PCPatchSetComputeOperator(), PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo() 643 @*/ 644 PetscErrorCode PCPatchSetComputeFunctionInteriorFacets(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 645 { 646 PC_PATCH *patch = (PC_PATCH *) pc->data; 647 648 PetscFunctionBegin; 649 patch->usercomputefintfacet = func; 650 patch->usercomputefintfacetctx = ctx; 651 PetscFunctionReturn(0); 652 } 653 654 /*@C 655 656 PCPatchSetComputeOperator - Set the callback used to compute patch matrices 657 658 Logically collective on PC 659 660 Input Parameters: 661 + pc - The PC 662 . func - The callback 663 - ctx - The user context 664 665 Calling sequence of func: 666 $ func (PC pc,PetscInt point,Vec x,Mat mat,IS facetIS,PetscInt n,const PetscInt* dofsArray,const PetscInt* dofsArrayWithAll,void* ctx) 667 668 + pc - The PC 669 . point - The point 670 . x - The input solution (not used in linear problems) 671 . mat - The patch matrix 672 . cellIS - An array of the cell numbers 673 . n - The size of dofsArray 674 . dofsArray - The dofmap for the dofs to be solved for 675 . dofsArrayWithAll - The dofmap for all dofs on the patch 676 - ctx - The user context 677 678 Level: advanced 679 680 Notes: 681 The matrix entries have been set to zero before the call. 682 683 .seealso: PCPatchGetComputeOperator(), PCPatchSetComputeFunction(), PCPatchSetDiscretisationInfo() 684 @*/ 685 PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 686 { 687 PC_PATCH *patch = (PC_PATCH *) pc->data; 688 689 PetscFunctionBegin; 690 patch->usercomputeop = func; 691 patch->usercomputeopctx = ctx; 692 PetscFunctionReturn(0); 693 } 694 695 /*@C 696 697 PCPatchSetComputeOperatorInteriorFacets - Set the callback used to compute facet integrals for patch matrices 698 699 Input Parameters: 700 + pc - The PC 701 . func - The callback 702 - ctx - The user context 703 704 Calling sequence of func: 705 $ func (PC pc,PetscInt point,Vec x,Mat mat,IS facetIS,PetscInt n,const PetscInt* dofsArray,const PetscInt* dofsArrayWithAll,void* ctx) 706 707 + pc - The PC 708 . point - The point 709 . x - The input solution (not used in linear problems) 710 . mat - The patch matrix 711 . facetIS - An array of the facet numbers 712 . n - The size of dofsArray 713 . dofsArray - The dofmap for the dofs to be solved for 714 . dofsArrayWithAll - The dofmap for all dofs on the patch 715 - ctx - The user context 716 717 Level: advanced 718 719 Notes: 720 The matrix entries have been set to zero before the call. 721 722 .seealso: PCPatchGetComputeOperator(), PCPatchSetComputeFunction(), PCPatchSetDiscretisationInfo() 723 @*/ 724 PetscErrorCode PCPatchSetComputeOperatorInteriorFacets(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 725 { 726 PC_PATCH *patch = (PC_PATCH *) pc->data; 727 728 PetscFunctionBegin; 729 patch->usercomputeopintfacet = func; 730 patch->usercomputeopintfacetctx = ctx; 731 PetscFunctionReturn(0); 732 } 733 734 /* On entry, ht contains the topological entities whose dofs we are responsible for solving for; 735 on exit, cht contains all the topological entities we need to compute their residuals. 736 In full generality this should incorporate knowledge of the sparsity pattern of the matrix; 737 here we assume a standard FE sparsity pattern.*/ 738 /* TODO: Use DMPlexGetAdjacency() */ 739 static PetscErrorCode PCPatchCompleteCellPatch(PC pc, PetscHSetI ht, PetscHSetI cht) 740 { 741 DM dm; 742 PC_PATCH *patch = (PC_PATCH *) pc->data; 743 PetscHashIter hi; 744 PetscInt point; 745 PetscInt *star = NULL, *closure = NULL; 746 PetscInt ignoredim, iStart = 0, iEnd = -1, starSize, closureSize, si, ci; 747 PetscInt *fStar = NULL, *fClosure = NULL; 748 PetscInt fBegin, fEnd, fsi, fci, fStarSize, fClosureSize; 749 PetscErrorCode ierr; 750 751 PetscFunctionBegin; 752 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 753 ierr = DMPlexGetHeightStratum(dm, 1, &fBegin, &fEnd);CHKERRQ(ierr); 754 ierr = PCPatchGetIgnoreDim(pc, &ignoredim);CHKERRQ(ierr); 755 if (ignoredim >= 0) {ierr = DMPlexGetDepthStratum(dm, ignoredim, &iStart, &iEnd);CHKERRQ(ierr);} 756 ierr = PetscHSetIClear(cht);CHKERRQ(ierr); 757 PetscHashIterBegin(ht, hi); 758 while (!PetscHashIterAtEnd(ht, hi)) { 759 760 PetscHashIterGetKey(ht, hi, point); 761 PetscHashIterNext(ht, hi); 762 763 /* Loop over all the cells that this point connects to */ 764 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 765 for (si = 0; si < starSize*2; si += 2) { 766 const PetscInt ownedpoint = star[si]; 767 /* TODO Check for point in cht before running through closure again */ 768 /* now loop over all entities in the closure of that cell */ 769 ierr = DMPlexGetTransitiveClosure(dm, ownedpoint, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 770 for (ci = 0; ci < closureSize*2; ci += 2) { 771 const PetscInt seenpoint = closure[ci]; 772 if (ignoredim >= 0 && seenpoint >= iStart && seenpoint < iEnd) continue; 773 ierr = PetscHSetIAdd(cht, seenpoint);CHKERRQ(ierr); 774 /* Facet integrals couple dofs across facets, so in that case for each of 775 * the facets we need to add all dofs on the other side of the facet to 776 * the seen dofs. */ 777 if(patch->usercomputeopintfacet){ 778 if(fBegin <= seenpoint && seenpoint < fEnd){ 779 ierr = DMPlexGetTransitiveClosure(dm, seenpoint, PETSC_FALSE, &fStarSize, &fStar);CHKERRQ(ierr); 780 for (fsi = 0; fsi < fStarSize*2; fsi += 2) { 781 ierr = DMPlexGetTransitiveClosure(dm, fStar[fsi], PETSC_TRUE, &fClosureSize, &fClosure);CHKERRQ(ierr); 782 for (fci = 0; fci < fClosureSize*2; fci += 2) { 783 ierr = PetscHSetIAdd(cht, fClosure[fci]);CHKERRQ(ierr); 784 } 785 ierr = DMPlexRestoreTransitiveClosure(dm, fStar[fsi], PETSC_TRUE, NULL, &fClosure);CHKERRQ(ierr); 786 } 787 ierr = DMPlexRestoreTransitiveClosure(dm, seenpoint, PETSC_FALSE, NULL, &fStar);CHKERRQ(ierr); 788 } 789 } 790 } 791 ierr = DMPlexRestoreTransitiveClosure(dm, ownedpoint, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr); 792 } 793 ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, NULL, &star);CHKERRQ(ierr); 794 } 795 PetscFunctionReturn(0); 796 } 797 798 static PetscErrorCode PCPatchGetGlobalDofs(PC pc, PetscSection dofSection[], PetscInt f, PetscBool combined, PetscInt p, PetscInt *dof, PetscInt *off) 799 { 800 PetscErrorCode ierr; 801 802 PetscFunctionBegin; 803 if (combined) { 804 if (f < 0) { 805 if (dof) {ierr = PetscSectionGetDof(dofSection[0], p, dof);CHKERRQ(ierr);} 806 if (off) {ierr = PetscSectionGetOffset(dofSection[0], p, off);CHKERRQ(ierr);} 807 } else { 808 if (dof) {ierr = PetscSectionGetFieldDof(dofSection[0], p, f, dof);CHKERRQ(ierr);} 809 if (off) {ierr = PetscSectionGetFieldOffset(dofSection[0], p, f, off);CHKERRQ(ierr);} 810 } 811 } else { 812 if (f < 0) { 813 PC_PATCH *patch = (PC_PATCH *) pc->data; 814 PetscInt fdof, g; 815 816 if (dof) { 817 *dof = 0; 818 for (g = 0; g < patch->nsubspaces; ++g) { 819 ierr = PetscSectionGetDof(dofSection[g], p, &fdof);CHKERRQ(ierr); 820 *dof += fdof; 821 } 822 } 823 if (off) { 824 *off = 0; 825 for (g = 0; g < patch->nsubspaces; ++g) { 826 ierr = PetscSectionGetOffset(dofSection[g], p, &fdof);CHKERRQ(ierr); 827 *off += fdof; 828 } 829 } 830 } else { 831 if (dof) {ierr = PetscSectionGetDof(dofSection[f], p, dof);CHKERRQ(ierr);} 832 if (off) {ierr = PetscSectionGetOffset(dofSection[f], p, off);CHKERRQ(ierr);} 833 } 834 } 835 PetscFunctionReturn(0); 836 } 837 838 /* Given a hash table with a set of topological entities (pts), compute the degrees of 839 freedom in global concatenated numbering on those entities. 840 For Vanka smoothing, this needs to do something special: ignore dofs of the 841 constraint subspace on entities that aren't the base entity we're building the patch 842 around. */ 843 static PetscErrorCode PCPatchGetPointDofs(PC pc, PetscHSetI pts, PetscHSetI dofs, PetscInt base, PetscHSetI* subspaces_to_exclude) 844 { 845 PC_PATCH *patch = (PC_PATCH *) pc->data; 846 PetscHashIter hi; 847 PetscInt ldof, loff; 848 PetscInt k, p; 849 PetscErrorCode ierr; 850 851 PetscFunctionBegin; 852 ierr = PetscHSetIClear(dofs);CHKERRQ(ierr); 853 for (k = 0; k < patch->nsubspaces; ++k) { 854 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 855 PetscInt bs = patch->bs[k]; 856 PetscInt j, l; 857 858 if (subspaces_to_exclude != NULL) { 859 PetscBool should_exclude_k = PETSC_FALSE; 860 PetscHSetIHas(*subspaces_to_exclude, k, &should_exclude_k); 861 if (should_exclude_k) { 862 /* only get this subspace dofs at the base entity, not any others */ 863 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, base, &ldof, &loff);CHKERRQ(ierr); 864 if (0 == ldof) continue; 865 for (j = loff; j < ldof + loff; ++j) { 866 for (l = 0; l < bs; ++l) { 867 PetscInt dof = bs*j + l + subspaceOffset; 868 ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr); 869 } 870 } 871 continue; /* skip the other dofs of this subspace */ 872 } 873 } 874 875 PetscHashIterBegin(pts, hi); 876 while (!PetscHashIterAtEnd(pts, hi)) { 877 PetscHashIterGetKey(pts, hi, p); 878 PetscHashIterNext(pts, hi); 879 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, p, &ldof, &loff);CHKERRQ(ierr); 880 if (0 == ldof) continue; 881 for (j = loff; j < ldof + loff; ++j) { 882 for (l = 0; l < bs; ++l) { 883 PetscInt dof = bs*j + l + subspaceOffset; 884 ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr); 885 } 886 } 887 } 888 } 889 PetscFunctionReturn(0); 890 } 891 892 /* Given two hash tables A and B, compute the keys in B that are not in A, and put them in C */ 893 static PetscErrorCode PCPatchComputeSetDifference_Private(PetscHSetI A, PetscHSetI B, PetscHSetI C) 894 { 895 PetscHashIter hi; 896 PetscInt key; 897 PetscBool flg; 898 PetscErrorCode ierr; 899 900 PetscFunctionBegin; 901 ierr = PetscHSetIClear(C);CHKERRQ(ierr); 902 PetscHashIterBegin(B, hi); 903 while (!PetscHashIterAtEnd(B, hi)) { 904 PetscHashIterGetKey(B, hi, key); 905 PetscHashIterNext(B, hi); 906 ierr = PetscHSetIHas(A, key, &flg);CHKERRQ(ierr); 907 if (!flg) {ierr = PetscHSetIAdd(C, key);CHKERRQ(ierr);} 908 } 909 PetscFunctionReturn(0); 910 } 911 912 /* 913 * PCPatchCreateCellPatches - create patches. 914 * 915 * Input Parameters: 916 * + dm - The DMPlex object defining the mesh 917 * 918 * Output Parameters: 919 * + cellCounts - Section with counts of cells around each vertex 920 * . cells - IS of the cell point indices of cells in each patch 921 * . pointCounts - Section with counts of cells around each vertex 922 * - point - IS of the cell point indices of cells in each patch 923 */ 924 static PetscErrorCode PCPatchCreateCellPatches(PC pc) 925 { 926 PC_PATCH *patch = (PC_PATCH *) pc->data; 927 DMLabel ghost = NULL; 928 DM dm, plex; 929 PetscHSetI ht, cht; 930 PetscSection cellCounts, pointCounts, intFacetCounts, extFacetCounts; 931 PetscInt *cellsArray, *pointsArray, *intFacetsArray, *extFacetsArray, *intFacetsToPatchCell; 932 PetscInt numCells, numPoints, numIntFacets, numExtFacets; 933 const PetscInt *leaves; 934 PetscInt nleaves, pStart, pEnd, cStart, cEnd, vStart, vEnd, v; 935 PetscBool isFiredrake; 936 PetscErrorCode ierr; 937 938 PetscFunctionBegin; 939 /* Used to keep track of the cells in the patch. */ 940 ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 941 ierr = PetscHSetICreate(&cht);CHKERRQ(ierr); 942 943 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 944 if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch PC\n"); 945 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 946 ierr = DMPlexGetChart(plex, &pStart, &pEnd);CHKERRQ(ierr); 947 ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr); 948 949 if (patch->user_patches) { 950 ierr = patch->userpatchconstructionop(pc, &patch->npatch, &patch->userIS, &patch->iterationSet, patch->userpatchconstructctx);CHKERRQ(ierr); 951 vStart = 0; vEnd = patch->npatch; 952 } else if (patch->ctype == PC_PATCH_PARDECOMP) { 953 vStart = 0; vEnd = 1; 954 } else if (patch->codim < 0) { 955 if (patch->dim < 0) {ierr = DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd);CHKERRQ(ierr);} 956 else {ierr = DMPlexGetDepthStratum(plex, patch->dim, &vStart, &vEnd);CHKERRQ(ierr);} 957 } else {ierr = DMPlexGetHeightStratum(plex, patch->codim, &vStart, &vEnd);CHKERRQ(ierr);} 958 patch->npatch = vEnd - vStart; 959 960 /* These labels mark the owned points. We only create patches around points that this process owns. */ 961 ierr = DMHasLabel(dm, "pyop2_ghost", &isFiredrake);CHKERRQ(ierr); 962 if (isFiredrake) { 963 ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr); 964 ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr); 965 } else { 966 PetscSF sf; 967 968 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 969 ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 970 nleaves = PetscMax(nleaves, 0); 971 } 972 973 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->cellCounts);CHKERRQ(ierr); 974 ierr = PetscObjectSetName((PetscObject) patch->cellCounts, "Patch Cell Layout");CHKERRQ(ierr); 975 cellCounts = patch->cellCounts; 976 ierr = PetscSectionSetChart(cellCounts, vStart, vEnd);CHKERRQ(ierr); 977 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->pointCounts);CHKERRQ(ierr); 978 ierr = PetscObjectSetName((PetscObject) patch->pointCounts, "Patch Point Layout");CHKERRQ(ierr); 979 pointCounts = patch->pointCounts; 980 ierr = PetscSectionSetChart(pointCounts, vStart, vEnd);CHKERRQ(ierr); 981 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->extFacetCounts);CHKERRQ(ierr); 982 ierr = PetscObjectSetName((PetscObject) patch->extFacetCounts, "Patch Exterior Facet Layout");CHKERRQ(ierr); 983 extFacetCounts = patch->extFacetCounts; 984 ierr = PetscSectionSetChart(extFacetCounts, vStart, vEnd);CHKERRQ(ierr); 985 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->intFacetCounts);CHKERRQ(ierr); 986 ierr = PetscObjectSetName((PetscObject) patch->intFacetCounts, "Patch Interior Facet Layout");CHKERRQ(ierr); 987 intFacetCounts = patch->intFacetCounts; 988 ierr = PetscSectionSetChart(intFacetCounts, vStart, vEnd);CHKERRQ(ierr); 989 /* Count cells and points in the patch surrounding each entity */ 990 PetscInt fStart, fEnd; 991 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 992 for (v = vStart; v < vEnd; ++v) { 993 PetscHashIter hi; 994 PetscInt chtSize, loc = -1; 995 PetscBool flg; 996 997 if (!patch->user_patches && patch->ctype != PC_PATCH_PARDECOMP) { 998 if (ghost) {ierr = DMLabelHasPoint(ghost, v, &flg);CHKERRQ(ierr);} 999 else {ierr = PetscFindInt(v, nleaves, leaves, &loc);CHKERRQ(ierr); flg = loc >=0 ? PETSC_TRUE : PETSC_FALSE;} 1000 /* Not an owned entity, don't make a cell patch. */ 1001 if (flg) continue; 1002 } 1003 1004 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 1005 ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr); 1006 ierr = PetscHSetIGetSize(cht, &chtSize);CHKERRQ(ierr); 1007 /* empty patch, continue */ 1008 if (chtSize == 0) continue; 1009 1010 /* safe because size(cht) > 0 from above */ 1011 PetscHashIterBegin(cht, hi); 1012 while (!PetscHashIterAtEnd(cht, hi)) { 1013 PetscInt point, pdof; 1014 1015 PetscHashIterGetKey(cht, hi, point); 1016 if (fStart <= point && point < fEnd) { 1017 const PetscInt *support; 1018 PetscInt supportSize, p; 1019 PetscBool interior = PETSC_TRUE; 1020 ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 1021 ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr); 1022 if (supportSize == 1) { 1023 interior = PETSC_FALSE; 1024 } else { 1025 for (p = 0; p < supportSize; p++) { 1026 PetscBool found; 1027 /* FIXME: can I do this while iterating over cht? */ 1028 PetscHSetIHas(cht, support[p], &found); 1029 if (!found) { 1030 interior = PETSC_FALSE; 1031 break; 1032 } 1033 } 1034 } 1035 if (interior) { 1036 ierr = PetscSectionAddDof(intFacetCounts, v, 1);CHKERRQ(ierr); 1037 } else { 1038 ierr = PetscSectionAddDof(extFacetCounts, v, 1);CHKERRQ(ierr); 1039 } 1040 } 1041 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr); 1042 if (pdof) {ierr = PetscSectionAddDof(pointCounts, v, 1);CHKERRQ(ierr);} 1043 if (point >= cStart && point < cEnd) {ierr = PetscSectionAddDof(cellCounts, v, 1);CHKERRQ(ierr);} 1044 PetscHashIterNext(cht, hi); 1045 } 1046 } 1047 if (isFiredrake) {ierr = DMLabelDestroyIndex(ghost);CHKERRQ(ierr);} 1048 1049 ierr = PetscSectionSetUp(cellCounts);CHKERRQ(ierr); 1050 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 1051 ierr = PetscMalloc1(numCells, &cellsArray);CHKERRQ(ierr); 1052 ierr = PetscSectionSetUp(pointCounts);CHKERRQ(ierr); 1053 ierr = PetscSectionGetStorageSize(pointCounts, &numPoints);CHKERRQ(ierr); 1054 ierr = PetscMalloc1(numPoints, &pointsArray);CHKERRQ(ierr); 1055 1056 ierr = PetscSectionSetUp(intFacetCounts);CHKERRQ(ierr); 1057 ierr = PetscSectionSetUp(extFacetCounts);CHKERRQ(ierr); 1058 ierr = PetscSectionGetStorageSize(intFacetCounts, &numIntFacets);CHKERRQ(ierr); 1059 ierr = PetscSectionGetStorageSize(extFacetCounts, &numExtFacets);CHKERRQ(ierr); 1060 ierr = PetscMalloc1(numIntFacets, &intFacetsArray);CHKERRQ(ierr); 1061 ierr = PetscMalloc1(numIntFacets*2, &intFacetsToPatchCell);CHKERRQ(ierr); 1062 ierr = PetscMalloc1(numExtFacets, &extFacetsArray);CHKERRQ(ierr); 1063 1064 1065 /* Now that we know how much space we need, run through again and actually remember the cells. */ 1066 for (v = vStart; v < vEnd; v++ ) { 1067 PetscHashIter hi; 1068 PetscInt dof, off, cdof, coff, efdof, efoff, ifdof, ifoff, pdof, n = 0, cn = 0, ifn = 0, efn = 0; 1069 1070 ierr = PetscSectionGetDof(pointCounts, v, &dof);CHKERRQ(ierr); 1071 ierr = PetscSectionGetOffset(pointCounts, v, &off);CHKERRQ(ierr); 1072 ierr = PetscSectionGetDof(cellCounts, v, &cdof);CHKERRQ(ierr); 1073 ierr = PetscSectionGetOffset(cellCounts, v, &coff);CHKERRQ(ierr); 1074 ierr = PetscSectionGetDof(intFacetCounts, v, &ifdof);CHKERRQ(ierr); 1075 ierr = PetscSectionGetOffset(intFacetCounts, v, &ifoff);CHKERRQ(ierr); 1076 ierr = PetscSectionGetDof(extFacetCounts, v, &efdof);CHKERRQ(ierr); 1077 ierr = PetscSectionGetOffset(extFacetCounts, v, &efoff);CHKERRQ(ierr); 1078 if (dof <= 0) continue; 1079 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 1080 ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr); 1081 PetscHashIterBegin(cht, hi); 1082 while (!PetscHashIterAtEnd(cht, hi)) { 1083 PetscInt point; 1084 1085 PetscHashIterGetKey(cht, hi, point); 1086 if (fStart <= point && point < fEnd) { 1087 const PetscInt *support; 1088 PetscInt supportSize, p; 1089 PetscBool interior = PETSC_TRUE; 1090 ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 1091 ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr); 1092 if (supportSize == 1) { 1093 interior = PETSC_FALSE; 1094 } else { 1095 for (p = 0; p < supportSize; p++) { 1096 PetscBool found; 1097 /* FIXME: can I do this while iterating over cht? */ 1098 PetscHSetIHas(cht, support[p], &found); 1099 if (!found) { 1100 interior = PETSC_FALSE; 1101 break; 1102 } 1103 } 1104 } 1105 if (interior) { 1106 intFacetsToPatchCell[2*(ifoff + ifn)] = support[0]; 1107 intFacetsToPatchCell[2*(ifoff + ifn) + 1] = support[1]; 1108 intFacetsArray[ifoff + ifn++] = point; 1109 } else { 1110 extFacetsArray[efoff + efn++] = point; 1111 } 1112 } 1113 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr); 1114 if (pdof) {pointsArray[off + n++] = point;} 1115 if (point >= cStart && point < cEnd) {cellsArray[coff + cn++] = point;} 1116 PetscHashIterNext(cht, hi); 1117 } 1118 if (ifn != ifdof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of interior facets in patch %D is %D, but should be %D", v, ifn, ifdof); 1119 if (efn != efdof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of exterior facets in patch %D is %D, but should be %D", v, efn, efdof); 1120 if (cn != cdof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of cells in patch %D is %D, but should be %D", v, cn, cdof); 1121 if (n != dof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of points in patch %D is %D, but should be %D", v, n, dof); 1122 1123 for (ifn = 0; ifn < ifdof; ifn++) { 1124 PetscInt cell0 = intFacetsToPatchCell[2*(ifoff + ifn)]; 1125 PetscInt cell1 = intFacetsToPatchCell[2*(ifoff + ifn) + 1]; 1126 PetscBool found0 = PETSC_FALSE, found1 = PETSC_FALSE; 1127 for (n = 0; n < cdof; n++) { 1128 if (!found0 && cell0 == cellsArray[coff + n]) { 1129 intFacetsToPatchCell[2*(ifoff + ifn)] = n; 1130 found0 = PETSC_TRUE; 1131 } 1132 if (!found1 && cell1 == cellsArray[coff + n]) { 1133 intFacetsToPatchCell[2*(ifoff + ifn) + 1] = n; 1134 found1 = PETSC_TRUE; 1135 } 1136 if (found0 && found1) break; 1137 } 1138 if (!(found0 && found1)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Didn't manage to find local point numbers for facet support"); 1139 } 1140 } 1141 ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 1142 ierr = PetscHSetIDestroy(&cht);CHKERRQ(ierr); 1143 ierr = DMDestroy(&plex);CHKERRQ(ierr); 1144 1145 ierr = ISCreateGeneral(PETSC_COMM_SELF, numCells, cellsArray, PETSC_OWN_POINTER, &patch->cells);CHKERRQ(ierr); 1146 ierr = PetscObjectSetName((PetscObject) patch->cells, "Patch Cells");CHKERRQ(ierr); 1147 if (patch->viewCells) { 1148 ierr = ObjectView((PetscObject) patch->cellCounts, patch->viewerCells, patch->formatCells);CHKERRQ(ierr); 1149 ierr = ObjectView((PetscObject) patch->cells, patch->viewerCells, patch->formatCells);CHKERRQ(ierr); 1150 } 1151 ierr = ISCreateGeneral(PETSC_COMM_SELF, numIntFacets, intFacetsArray, PETSC_OWN_POINTER, &patch->intFacets);CHKERRQ(ierr); 1152 ierr = PetscObjectSetName((PetscObject) patch->intFacets, "Patch Interior Facets");CHKERRQ(ierr); 1153 ierr = ISCreateGeneral(PETSC_COMM_SELF, 2*numIntFacets, intFacetsToPatchCell, PETSC_OWN_POINTER, &patch->intFacetsToPatchCell);CHKERRQ(ierr); 1154 ierr = PetscObjectSetName((PetscObject) patch->intFacetsToPatchCell, "Patch Interior Facets local support");CHKERRQ(ierr); 1155 if (patch->viewIntFacets) { 1156 ierr = ObjectView((PetscObject) patch->intFacetCounts, patch->viewerIntFacets, patch->formatIntFacets);CHKERRQ(ierr); 1157 ierr = ObjectView((PetscObject) patch->intFacets, patch->viewerIntFacets, patch->formatIntFacets);CHKERRQ(ierr); 1158 ierr = ObjectView((PetscObject) patch->intFacetsToPatchCell, patch->viewerIntFacets, patch->formatIntFacets);CHKERRQ(ierr); 1159 } 1160 ierr = ISCreateGeneral(PETSC_COMM_SELF, numExtFacets, extFacetsArray, PETSC_OWN_POINTER, &patch->extFacets);CHKERRQ(ierr); 1161 ierr = PetscObjectSetName((PetscObject) patch->extFacets, "Patch Exterior Facets");CHKERRQ(ierr); 1162 if (patch->viewExtFacets) { 1163 ierr = ObjectView((PetscObject) patch->extFacetCounts, patch->viewerExtFacets, patch->formatExtFacets);CHKERRQ(ierr); 1164 ierr = ObjectView((PetscObject) patch->extFacets, patch->viewerExtFacets, patch->formatExtFacets);CHKERRQ(ierr); 1165 } 1166 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, pointsArray, PETSC_OWN_POINTER, &patch->points);CHKERRQ(ierr); 1167 ierr = PetscObjectSetName((PetscObject) patch->points, "Patch Points");CHKERRQ(ierr); 1168 if (patch->viewPoints) { 1169 ierr = ObjectView((PetscObject) patch->pointCounts, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr); 1170 ierr = ObjectView((PetscObject) patch->points, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr); 1171 } 1172 PetscFunctionReturn(0); 1173 } 1174 1175 /* 1176 * PCPatchCreateCellPatchDiscretisationInfo - Build the dof maps for cell patches 1177 * 1178 * Input Parameters: 1179 * + dm - The DMPlex object defining the mesh 1180 * . cellCounts - Section with counts of cells around each vertex 1181 * . cells - IS of the cell point indices of cells in each patch 1182 * . cellNumbering - Section mapping plex cell points to Firedrake cell indices. 1183 * . nodesPerCell - number of nodes per cell. 1184 * - cellNodeMap - map from cells to node indices (nodesPerCell * numCells) 1185 * 1186 * Output Parameters: 1187 * + dofs - IS of local dof numbers of each cell in the patch, where local is a patch local numbering 1188 * . gtolCounts - Section with counts of dofs per cell patch 1189 * - gtol - IS mapping from global dofs to local dofs for each patch. 1190 */ 1191 static PetscErrorCode PCPatchCreateCellPatchDiscretisationInfo(PC pc) 1192 { 1193 PC_PATCH *patch = (PC_PATCH *) pc->data; 1194 PetscSection cellCounts = patch->cellCounts; 1195 PetscSection pointCounts = patch->pointCounts; 1196 PetscSection gtolCounts, gtolCountsWithArtificial = NULL, gtolCountsWithAll = NULL; 1197 IS cells = patch->cells; 1198 IS points = patch->points; 1199 PetscSection cellNumbering = patch->cellNumbering; 1200 PetscInt Nf = patch->nsubspaces; 1201 PetscInt numCells, numPoints; 1202 PetscInt numDofs; 1203 PetscInt numGlobalDofs, numGlobalDofsWithArtificial, numGlobalDofsWithAll; 1204 PetscInt totalDofsPerCell = patch->totalDofsPerCell; 1205 PetscInt vStart, vEnd, v; 1206 const PetscInt *cellsArray, *pointsArray; 1207 PetscInt *newCellsArray = NULL; 1208 PetscInt *dofsArray = NULL; 1209 PetscInt *dofsArrayWithArtificial = NULL; 1210 PetscInt *dofsArrayWithAll = NULL; 1211 PetscInt *offsArray = NULL; 1212 PetscInt *offsArrayWithArtificial = NULL; 1213 PetscInt *offsArrayWithAll = NULL; 1214 PetscInt *asmArray = NULL; 1215 PetscInt *asmArrayWithArtificial = NULL; 1216 PetscInt *asmArrayWithAll = NULL; 1217 PetscInt *globalDofsArray = NULL; 1218 PetscInt *globalDofsArrayWithArtificial = NULL; 1219 PetscInt *globalDofsArrayWithAll = NULL; 1220 PetscInt globalIndex = 0; 1221 PetscInt key = 0; 1222 PetscInt asmKey = 0; 1223 DM dm = NULL; 1224 const PetscInt *bcNodes = NULL; 1225 PetscHMapI ht; 1226 PetscHMapI htWithArtificial; 1227 PetscHMapI htWithAll; 1228 PetscHSetI globalBcs; 1229 PetscInt numBcs; 1230 PetscHSetI ownedpts, seenpts, owneddofs, seendofs, artificialbcs; 1231 PetscInt pStart, pEnd, p, i; 1232 char option[PETSC_MAX_PATH_LEN]; 1233 PetscBool isNonlinear; 1234 PetscErrorCode ierr; 1235 1236 PetscFunctionBegin; 1237 1238 ierr = PCGetDM(pc, &dm); CHKERRQ(ierr); 1239 /* dofcounts section is cellcounts section * dofPerCell */ 1240 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 1241 ierr = PetscSectionGetStorageSize(patch->pointCounts, &numPoints);CHKERRQ(ierr); 1242 numDofs = numCells * totalDofsPerCell; 1243 ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr); 1244 ierr = PetscMalloc1(numPoints*Nf, &offsArray);CHKERRQ(ierr); 1245 ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr); 1246 ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr); 1247 ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr); 1248 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr); 1249 gtolCounts = patch->gtolCounts; 1250 ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr); 1251 ierr = PetscObjectSetName((PetscObject) patch->gtolCounts, "Patch Global Index Section");CHKERRQ(ierr); 1252 1253 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) 1254 { 1255 ierr = PetscMalloc1(numPoints*Nf, &offsArrayWithArtificial);CHKERRQ(ierr); 1256 ierr = PetscMalloc1(numDofs, &asmArrayWithArtificial);CHKERRQ(ierr); 1257 ierr = PetscMalloc1(numDofs, &dofsArrayWithArtificial);CHKERRQ(ierr); 1258 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCountsWithArtificial);CHKERRQ(ierr); 1259 gtolCountsWithArtificial = patch->gtolCountsWithArtificial; 1260 ierr = PetscSectionSetChart(gtolCountsWithArtificial, vStart, vEnd);CHKERRQ(ierr); 1261 ierr = PetscObjectSetName((PetscObject) patch->gtolCountsWithArtificial, "Patch Global Index Section Including Artificial BCs");CHKERRQ(ierr); 1262 } 1263 1264 isNonlinear = patch->isNonlinear; 1265 if(isNonlinear) 1266 { 1267 ierr = PetscMalloc1(numPoints*Nf, &offsArrayWithAll);CHKERRQ(ierr); 1268 ierr = PetscMalloc1(numDofs, &asmArrayWithAll);CHKERRQ(ierr); 1269 ierr = PetscMalloc1(numDofs, &dofsArrayWithAll);CHKERRQ(ierr); 1270 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCountsWithAll);CHKERRQ(ierr); 1271 gtolCountsWithAll = patch->gtolCountsWithAll; 1272 ierr = PetscSectionSetChart(gtolCountsWithAll, vStart, vEnd);CHKERRQ(ierr); 1273 ierr = PetscObjectSetName((PetscObject) patch->gtolCountsWithAll, "Patch Global Index Section Including All BCs");CHKERRQ(ierr); 1274 } 1275 1276 /* Outside the patch loop, get the dofs that are globally-enforced Dirichlet 1277 conditions */ 1278 ierr = PetscHSetICreate(&globalBcs);CHKERRQ(ierr); 1279 ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr); 1280 ierr = ISGetSize(patch->ghostBcNodes, &numBcs); CHKERRQ(ierr); 1281 for (i = 0; i < numBcs; ++i) { 1282 ierr = PetscHSetIAdd(globalBcs, bcNodes[i]);CHKERRQ(ierr); /* these are already in concatenated numbering */ 1283 } 1284 ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr); 1285 ierr = ISDestroy(&patch->ghostBcNodes); CHKERRQ(ierr); /* memory optimisation */ 1286 1287 /* Hash tables for artificial BC construction */ 1288 ierr = PetscHSetICreate(&ownedpts);CHKERRQ(ierr); 1289 ierr = PetscHSetICreate(&seenpts);CHKERRQ(ierr); 1290 ierr = PetscHSetICreate(&owneddofs);CHKERRQ(ierr); 1291 ierr = PetscHSetICreate(&seendofs);CHKERRQ(ierr); 1292 ierr = PetscHSetICreate(&artificialbcs);CHKERRQ(ierr); 1293 1294 ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr); 1295 ierr = ISGetIndices(points, &pointsArray);CHKERRQ(ierr); 1296 ierr = PetscHMapICreate(&ht);CHKERRQ(ierr); 1297 ierr = PetscHMapICreate(&htWithArtificial);CHKERRQ(ierr); 1298 ierr = PetscHMapICreate(&htWithAll);CHKERRQ(ierr); 1299 for (v = vStart; v < vEnd; ++v) { 1300 PetscInt localIndex = 0; 1301 PetscInt localIndexWithArtificial = 0; 1302 PetscInt localIndexWithAll = 0; 1303 PetscInt dof, off, i, j, k, l; 1304 1305 ierr = PetscHMapIClear(ht);CHKERRQ(ierr); 1306 ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr); 1307 ierr = PetscHMapIClear(htWithAll);CHKERRQ(ierr); 1308 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 1309 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 1310 if (dof <= 0) continue; 1311 1312 /* Calculate the global numbers of the artificial BC dofs here first */ 1313 ierr = patch->patchconstructop((void*)patch, dm, v, ownedpts); CHKERRQ(ierr); 1314 ierr = PCPatchCompleteCellPatch(pc, ownedpts, seenpts); CHKERRQ(ierr); 1315 ierr = PCPatchGetPointDofs(pc, ownedpts, owneddofs, v, &patch->subspaces_to_exclude); CHKERRQ(ierr); 1316 ierr = PCPatchGetPointDofs(pc, seenpts, seendofs, v, NULL); CHKERRQ(ierr); 1317 ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs); CHKERRQ(ierr); 1318 if (patch->viewPatches) { 1319 PetscHSetI globalbcdofs; 1320 PetscHashIter hi; 1321 MPI_Comm comm = PetscObjectComm((PetscObject)pc); 1322 1323 ierr = PetscHSetICreate(&globalbcdofs);CHKERRQ(ierr); 1324 ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v); CHKERRQ(ierr); 1325 PetscHashIterBegin(owneddofs, hi); 1326 while (!PetscHashIterAtEnd(owneddofs, hi)) { 1327 PetscInt globalDof; 1328 1329 PetscHashIterGetKey(owneddofs, hi, globalDof); 1330 PetscHashIterNext(owneddofs, hi); 1331 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 1332 } 1333 ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr); 1334 ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v); CHKERRQ(ierr); 1335 PetscHashIterBegin(seendofs, hi); 1336 while (!PetscHashIterAtEnd(seendofs, hi)) { 1337 PetscInt globalDof; 1338 PetscBool flg; 1339 1340 PetscHashIterGetKey(seendofs, hi, globalDof); 1341 PetscHashIterNext(seendofs, hi); 1342 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 1343 1344 ierr = PetscHSetIHas(globalBcs, globalDof, &flg);CHKERRQ(ierr); 1345 if (flg) {ierr = PetscHSetIAdd(globalbcdofs, globalDof);CHKERRQ(ierr);} 1346 } 1347 ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr); 1348 ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v); CHKERRQ(ierr); 1349 ierr = PetscHSetIGetSize(globalbcdofs, &numBcs);CHKERRQ(ierr); 1350 if (numBcs > 0) { 1351 PetscHashIterBegin(globalbcdofs, hi); 1352 while (!PetscHashIterAtEnd(globalbcdofs, hi)) { 1353 PetscInt globalDof; 1354 PetscHashIterGetKey(globalbcdofs, hi, globalDof); 1355 PetscHashIterNext(globalbcdofs, hi); 1356 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr); 1357 } 1358 } 1359 ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr); 1360 ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr); 1361 ierr = PetscHSetIGetSize(artificialbcs, &numBcs);CHKERRQ(ierr); 1362 if (numBcs > 0) { 1363 PetscHashIterBegin(artificialbcs, hi); 1364 while (!PetscHashIterAtEnd(artificialbcs, hi)) { 1365 PetscInt globalDof; 1366 PetscHashIterGetKey(artificialbcs, hi, globalDof); 1367 PetscHashIterNext(artificialbcs, hi); 1368 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 1369 } 1370 } 1371 ierr = PetscSynchronizedPrintf(comm, "\n\n"); CHKERRQ(ierr); 1372 ierr = PetscHSetIDestroy(&globalbcdofs);CHKERRQ(ierr); 1373 } 1374 for (k = 0; k < patch->nsubspaces; ++k) { 1375 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1376 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1377 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1378 PetscInt bs = patch->bs[k]; 1379 1380 for (i = off; i < off + dof; ++i) { 1381 /* Walk over the cells in this patch. */ 1382 const PetscInt c = cellsArray[i]; 1383 PetscInt cell = c; 1384 1385 /* TODO Change this to an IS */ 1386 if (cellNumbering) { 1387 ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr); 1388 if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c); 1389 ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr); 1390 } 1391 newCellsArray[i] = cell; 1392 for (j = 0; j < nodesPerCell; ++j) { 1393 /* For each global dof, map it into contiguous local storage. */ 1394 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset; 1395 /* finally, loop over block size */ 1396 for (l = 0; l < bs; ++l) { 1397 PetscInt localDof; 1398 PetscBool isGlobalBcDof, isArtificialBcDof; 1399 1400 /* first, check if this is either a globally enforced or locally enforced BC dof */ 1401 ierr = PetscHSetIHas(globalBcs, globalDof + l, &isGlobalBcDof);CHKERRQ(ierr); 1402 ierr = PetscHSetIHas(artificialbcs, globalDof + l, &isArtificialBcDof);CHKERRQ(ierr); 1403 1404 /* if it's either, don't ever give it a local dof number */ 1405 if (isGlobalBcDof || isArtificialBcDof) { 1406 dofsArray[globalIndex] = -1; /* don't use this in assembly in this patch */ 1407 } else { 1408 ierr = PetscHMapIGet(ht, globalDof + l, &localDof);CHKERRQ(ierr); 1409 if (localDof == -1) { 1410 localDof = localIndex++; 1411 ierr = PetscHMapISet(ht, globalDof + l, localDof);CHKERRQ(ierr); 1412 } 1413 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 1414 /* And store. */ 1415 dofsArray[globalIndex] = localDof; 1416 } 1417 1418 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1419 if (isGlobalBcDof) { 1420 dofsArrayWithArtificial[globalIndex] = -1; /* don't use this in assembly in this patch */ 1421 } else { 1422 ierr = PetscHMapIGet(htWithArtificial, globalDof + l, &localDof);CHKERRQ(ierr); 1423 if (localDof == -1) { 1424 localDof = localIndexWithArtificial++; 1425 ierr = PetscHMapISet(htWithArtificial, globalDof + l, localDof);CHKERRQ(ierr); 1426 } 1427 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 1428 /* And store.*/ 1429 dofsArrayWithArtificial[globalIndex] = localDof; 1430 } 1431 } 1432 1433 if(isNonlinear) { 1434 /* Build the dofmap for the function space with _all_ dofs, 1435 including those in any kind of boundary condition */ 1436 ierr = PetscHMapIGet(htWithAll, globalDof + l, &localDof);CHKERRQ(ierr); 1437 if (localDof == -1) { 1438 localDof = localIndexWithAll++; 1439 ierr = PetscHMapISet(htWithAll, globalDof + l, localDof);CHKERRQ(ierr); 1440 } 1441 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 1442 /* And store.*/ 1443 dofsArrayWithAll[globalIndex] = localDof; 1444 } 1445 globalIndex++; 1446 } 1447 } 1448 } 1449 } 1450 /*How many local dofs in this patch? */ 1451 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1452 ierr = PetscHMapIGetSize(htWithArtificial, &dof);CHKERRQ(ierr); 1453 ierr = PetscSectionSetDof(gtolCountsWithArtificial, v, dof);CHKERRQ(ierr); 1454 } 1455 if (isNonlinear) { 1456 ierr = PetscHMapIGetSize(htWithAll, &dof);CHKERRQ(ierr); 1457 ierr = PetscSectionSetDof(gtolCountsWithAll, v, dof);CHKERRQ(ierr); 1458 } 1459 ierr = PetscHMapIGetSize(ht, &dof);CHKERRQ(ierr); 1460 ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr); 1461 } 1462 if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex); 1463 ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr); 1464 ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr); 1465 ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr); 1466 1467 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1468 ierr = PetscSectionSetUp(gtolCountsWithArtificial);CHKERRQ(ierr); 1469 ierr = PetscSectionGetStorageSize(gtolCountsWithArtificial, &numGlobalDofsWithArtificial);CHKERRQ(ierr); 1470 ierr = PetscMalloc1(numGlobalDofsWithArtificial, &globalDofsArrayWithArtificial);CHKERRQ(ierr); 1471 } 1472 if (isNonlinear) { 1473 ierr = PetscSectionSetUp(gtolCountsWithAll);CHKERRQ(ierr); 1474 ierr = PetscSectionGetStorageSize(gtolCountsWithAll, &numGlobalDofsWithAll);CHKERRQ(ierr); 1475 ierr = PetscMalloc1(numGlobalDofsWithAll, &globalDofsArrayWithAll);CHKERRQ(ierr); 1476 } 1477 /* Now populate the global to local map. This could be merged into the above loop if we were willing to deal with reallocs. */ 1478 for (v = vStart; v < vEnd; ++v) { 1479 PetscHashIter hi; 1480 PetscInt dof, off, Np, ooff, i, j, k, l; 1481 1482 ierr = PetscHMapIClear(ht);CHKERRQ(ierr); 1483 ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr); 1484 ierr = PetscHMapIClear(htWithAll);CHKERRQ(ierr); 1485 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 1486 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 1487 ierr = PetscSectionGetDof(pointCounts, v, &Np);CHKERRQ(ierr); 1488 ierr = PetscSectionGetOffset(pointCounts, v, &ooff);CHKERRQ(ierr); 1489 if (dof <= 0) continue; 1490 1491 for (k = 0; k < patch->nsubspaces; ++k) { 1492 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1493 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1494 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1495 PetscInt bs = patch->bs[k]; 1496 PetscInt goff; 1497 1498 for (i = off; i < off + dof; ++i) { 1499 /* Reconstruct mapping of global-to-local on this patch. */ 1500 const PetscInt c = cellsArray[i]; 1501 PetscInt cell = c; 1502 1503 if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);} 1504 for (j = 0; j < nodesPerCell; ++j) { 1505 for (l = 0; l < bs; ++l) { 1506 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset; 1507 const PetscInt localDof = dofsArray[key]; 1508 if (localDof >= 0) {ierr = PetscHMapISet(ht, globalDof, localDof);CHKERRQ(ierr);} 1509 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1510 const PetscInt localDofWithArtificial = dofsArrayWithArtificial[key]; 1511 if (localDofWithArtificial >= 0) { 1512 ierr = PetscHMapISet(htWithArtificial, globalDof, localDofWithArtificial);CHKERRQ(ierr); 1513 } 1514 } 1515 if (isNonlinear) { 1516 const PetscInt localDofWithAll = dofsArrayWithAll[key]; 1517 if (localDofWithAll >= 0) { 1518 ierr = PetscHMapISet(htWithAll, globalDof, localDofWithAll);CHKERRQ(ierr); 1519 } 1520 } 1521 key++; 1522 } 1523 } 1524 } 1525 1526 /* Shove it in the output data structure. */ 1527 ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr); 1528 PetscHashIterBegin(ht, hi); 1529 while (!PetscHashIterAtEnd(ht, hi)) { 1530 PetscInt globalDof, localDof; 1531 1532 PetscHashIterGetKey(ht, hi, globalDof); 1533 PetscHashIterGetVal(ht, hi, localDof); 1534 if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof; 1535 PetscHashIterNext(ht, hi); 1536 } 1537 1538 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1539 ierr = PetscSectionGetOffset(gtolCountsWithArtificial, v, &goff);CHKERRQ(ierr); 1540 PetscHashIterBegin(htWithArtificial, hi); 1541 while (!PetscHashIterAtEnd(htWithArtificial, hi)) { 1542 PetscInt globalDof, localDof; 1543 PetscHashIterGetKey(htWithArtificial, hi, globalDof); 1544 PetscHashIterGetVal(htWithArtificial, hi, localDof); 1545 if (globalDof >= 0) globalDofsArrayWithArtificial[goff + localDof] = globalDof; 1546 PetscHashIterNext(htWithArtificial, hi); 1547 } 1548 } 1549 if (isNonlinear) { 1550 ierr = PetscSectionGetOffset(gtolCountsWithAll, v, &goff);CHKERRQ(ierr); 1551 PetscHashIterBegin(htWithAll, hi); 1552 while (!PetscHashIterAtEnd(htWithAll, hi)) { 1553 PetscInt globalDof, localDof; 1554 PetscHashIterGetKey(htWithAll, hi, globalDof); 1555 PetscHashIterGetVal(htWithAll, hi, localDof); 1556 if (globalDof >= 0) globalDofsArrayWithAll[goff + localDof] = globalDof; 1557 PetscHashIterNext(htWithAll, hi); 1558 } 1559 } 1560 1561 for (p = 0; p < Np; ++p) { 1562 const PetscInt point = pointsArray[ooff + p]; 1563 PetscInt globalDof, localDof; 1564 1565 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, point, NULL, &globalDof);CHKERRQ(ierr); 1566 ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr); 1567 offsArray[(ooff + p)*Nf + k] = localDof; 1568 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1569 ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr); 1570 offsArrayWithArtificial[(ooff + p)*Nf + k] = localDof; 1571 } 1572 if (isNonlinear) { 1573 ierr = PetscHMapIGet(htWithAll, globalDof, &localDof);CHKERRQ(ierr); 1574 offsArrayWithAll[(ooff + p)*Nf + k] = localDof; 1575 } 1576 } 1577 } 1578 1579 ierr = PetscHSetIDestroy(&globalBcs);CHKERRQ(ierr); 1580 ierr = PetscHSetIDestroy(&ownedpts);CHKERRQ(ierr); 1581 ierr = PetscHSetIDestroy(&seenpts);CHKERRQ(ierr); 1582 ierr = PetscHSetIDestroy(&owneddofs);CHKERRQ(ierr); 1583 ierr = PetscHSetIDestroy(&seendofs);CHKERRQ(ierr); 1584 ierr = PetscHSetIDestroy(&artificialbcs);CHKERRQ(ierr); 1585 1586 /* At this point, we have a hash table ht built that maps globalDof -> localDof. 1587 We need to create the dof table laid out cellwise first, then by subspace, 1588 as the assembler assembles cell-wise and we need to stuff the different 1589 contributions of the different function spaces to the right places. So we loop 1590 over cells, then over subspaces. */ 1591 if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */ 1592 for (i = off; i < off + dof; ++i) { 1593 const PetscInt c = cellsArray[i]; 1594 PetscInt cell = c; 1595 1596 if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);} 1597 for (k = 0; k < patch->nsubspaces; ++k) { 1598 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1599 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1600 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1601 PetscInt bs = patch->bs[k]; 1602 1603 for (j = 0; j < nodesPerCell; ++j) { 1604 for (l = 0; l < bs; ++l) { 1605 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset; 1606 PetscInt localDof; 1607 1608 ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr); 1609 /* If it's not in the hash table, i.e. is a BC dof, 1610 then the PetscHSetIMap above gives -1, which matches 1611 exactly the convention for PETSc's matrix assembly to 1612 ignore the dof. So we don't need to do anything here */ 1613 asmArray[asmKey] = localDof; 1614 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1615 ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr); 1616 asmArrayWithArtificial[asmKey] = localDof; 1617 } 1618 if (isNonlinear) { 1619 ierr = PetscHMapIGet(htWithAll, globalDof, &localDof);CHKERRQ(ierr); 1620 asmArrayWithAll[asmKey] = localDof; 1621 } 1622 asmKey++; 1623 } 1624 } 1625 } 1626 } 1627 } 1628 } 1629 if (1 == patch->nsubspaces) { 1630 ierr = PetscMemcpy(asmArray, dofsArray, numDofs * sizeof(PetscInt));CHKERRQ(ierr); 1631 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1632 ierr = PetscMemcpy(asmArrayWithArtificial, dofsArrayWithArtificial, numDofs * sizeof(PetscInt));CHKERRQ(ierr); 1633 } 1634 if (isNonlinear) { 1635 ierr = PetscMemcpy(asmArrayWithAll, dofsArrayWithAll, numDofs * sizeof(PetscInt));CHKERRQ(ierr); 1636 } 1637 } 1638 1639 ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr); 1640 ierr = PetscHMapIDestroy(&htWithArtificial);CHKERRQ(ierr); 1641 ierr = PetscHMapIDestroy(&htWithAll);CHKERRQ(ierr); 1642 ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr); 1643 ierr = ISRestoreIndices(points, &pointsArray);CHKERRQ(ierr); 1644 ierr = PetscFree(dofsArray);CHKERRQ(ierr); 1645 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1646 ierr = PetscFree(dofsArrayWithArtificial);CHKERRQ(ierr); 1647 } 1648 if (isNonlinear) { 1649 ierr = PetscFree(dofsArrayWithAll);CHKERRQ(ierr); 1650 } 1651 /* Create placeholder section for map from points to patch dofs */ 1652 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->patchSection);CHKERRQ(ierr); 1653 ierr = PetscSectionSetNumFields(patch->patchSection, patch->nsubspaces);CHKERRQ(ierr); 1654 if (patch->combined) { 1655 PetscInt numFields; 1656 ierr = PetscSectionGetNumFields(patch->dofSection[0], &numFields);CHKERRQ(ierr); 1657 if (numFields != patch->nsubspaces) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Mismatch between number of section fields %D and number of subspaces %D", numFields, patch->nsubspaces); 1658 ierr = PetscSectionGetChart(patch->dofSection[0], &pStart, &pEnd);CHKERRQ(ierr); 1659 ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr); 1660 for (p = pStart; p < pEnd; ++p) { 1661 PetscInt dof, fdof, f; 1662 1663 ierr = PetscSectionGetDof(patch->dofSection[0], p, &dof);CHKERRQ(ierr); 1664 ierr = PetscSectionSetDof(patch->patchSection, p, dof);CHKERRQ(ierr); 1665 for (f = 0; f < patch->nsubspaces; ++f) { 1666 ierr = PetscSectionGetFieldDof(patch->dofSection[0], p, f, &fdof);CHKERRQ(ierr); 1667 ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr); 1668 } 1669 } 1670 } else { 1671 PetscInt pStartf, pEndf, f; 1672 pStart = PETSC_MAX_INT; 1673 pEnd = PETSC_MIN_INT; 1674 for (f = 0; f < patch->nsubspaces; ++f) { 1675 ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr); 1676 pStart = PetscMin(pStart, pStartf); 1677 pEnd = PetscMax(pEnd, pEndf); 1678 } 1679 ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr); 1680 for (f = 0; f < patch->nsubspaces; ++f) { 1681 ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr); 1682 for (p = pStartf; p < pEndf; ++p) { 1683 PetscInt fdof; 1684 ierr = PetscSectionGetDof(patch->dofSection[f], p, &fdof);CHKERRQ(ierr); 1685 ierr = PetscSectionAddDof(patch->patchSection, p, fdof);CHKERRQ(ierr); 1686 ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr); 1687 } 1688 } 1689 } 1690 ierr = PetscSectionSetUp(patch->patchSection);CHKERRQ(ierr); 1691 ierr = PetscSectionSetUseFieldOffsets(patch->patchSection, PETSC_TRUE);CHKERRQ(ierr); 1692 /* Replace cell indices with firedrake-numbered ones. */ 1693 ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr); 1694 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr); 1695 ierr = PetscObjectSetName((PetscObject) patch->gtol, "Global Indices");CHKERRQ(ierr); 1696 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_g2l_view", patch->classname);CHKERRQ(ierr); 1697 ierr = PetscSectionViewFromOptions(patch->gtolCounts, (PetscObject) pc, option);CHKERRQ(ierr); 1698 ierr = ISViewFromOptions(patch->gtol, (PetscObject) pc, option);CHKERRQ(ierr); 1699 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr); 1700 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArray, PETSC_OWN_POINTER, &patch->offs);CHKERRQ(ierr); 1701 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1702 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofsWithArtificial, globalDofsArrayWithArtificial, PETSC_OWN_POINTER, &patch->gtolWithArtificial);CHKERRQ(ierr); 1703 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArrayWithArtificial, PETSC_OWN_POINTER, &patch->dofsWithArtificial);CHKERRQ(ierr); 1704 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArrayWithArtificial, PETSC_OWN_POINTER, &patch->offsWithArtificial);CHKERRQ(ierr); 1705 } 1706 if (isNonlinear) { 1707 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofsWithAll, globalDofsArrayWithAll, PETSC_OWN_POINTER, &patch->gtolWithAll);CHKERRQ(ierr); 1708 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArrayWithAll, PETSC_OWN_POINTER, &patch->dofsWithAll);CHKERRQ(ierr); 1709 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArrayWithAll, PETSC_OWN_POINTER, &patch->offsWithAll);CHKERRQ(ierr); 1710 } 1711 PetscFunctionReturn(0); 1712 } 1713 1714 static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat, PetscBool withArtificial) 1715 { 1716 PC_PATCH *patch = (PC_PATCH *) pc->data; 1717 Vec x, y; 1718 PetscBool flg; 1719 PetscInt csize, rsize; 1720 const char *prefix = NULL; 1721 PetscErrorCode ierr; 1722 1723 PetscFunctionBegin; 1724 if(withArtificial) { 1725 /* would be nice if we could create a rectangular matrix of size numDofsWithArtificial x numDofs here */ 1726 x = patch->patchRHSWithArtificial[point]; 1727 y = patch->patchRHSWithArtificial[point]; 1728 } else { 1729 x = patch->patchRHS[point]; 1730 y = patch->patchUpdate[point]; 1731 } 1732 1733 ierr = VecGetSize(x, &csize);CHKERRQ(ierr); 1734 ierr = VecGetSize(y, &rsize);CHKERRQ(ierr); 1735 ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr); 1736 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 1737 ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr); 1738 ierr = MatAppendOptionsPrefix(*mat, "pc_patch_sub_");CHKERRQ(ierr); 1739 if (patch->sub_mat_type) {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);} 1740 else if (!patch->sub_mat_type) {ierr = MatSetType(*mat, MATDENSE);CHKERRQ(ierr);} 1741 ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr); 1742 ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr); 1743 if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);} 1744 /* Sparse patch matrices */ 1745 if (!flg) { 1746 PetscBT bt; 1747 PetscInt *dnnz = NULL; 1748 const PetscInt *dofsArray = NULL; 1749 PetscInt pStart, pEnd, ncell, offset, c, i, j; 1750 1751 if(withArtificial) { 1752 ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 1753 } else { 1754 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1755 } 1756 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1757 point += pStart; 1758 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1759 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1760 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1761 ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1762 /* A PetscBT uses N^2 bits to store the sparsity pattern on a 1763 * patch. This is probably OK if the patches are not too big, 1764 * but uses too much memory. We therefore switch based on rsize. */ 1765 if (rsize < 3000) { /* FIXME: I picked this switch value out of my hat */ 1766 PetscScalar *zeroes; 1767 PetscInt rows; 1768 1769 ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr); 1770 ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr); 1771 for (c = 0; c < ncell; ++c) { 1772 const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell; 1773 for (i = 0; i < patch->totalDofsPerCell; ++i) { 1774 const PetscInt row = idx[i]; 1775 if (row < 0) continue; 1776 for (j = 0; j < patch->totalDofsPerCell; ++j) { 1777 const PetscInt col = idx[j]; 1778 const PetscInt key = row*rsize + col; 1779 if (col < 0) continue; 1780 if (!PetscBTLookupSet(bt, key)) ++dnnz[row]; 1781 } 1782 } 1783 } 1784 1785 if (patch->usercomputeopintfacet) { 1786 const PetscInt *intFacetsArray = NULL; 1787 PetscInt i, numIntFacets, intFacetOffset; 1788 const PetscInt *facetCells = NULL; 1789 1790 ierr = PetscSectionGetDof(patch->intFacetCounts, point, &numIntFacets);CHKERRQ(ierr); 1791 ierr = PetscSectionGetOffset(patch->intFacetCounts, point, &intFacetOffset);CHKERRQ(ierr); 1792 ierr = ISGetIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 1793 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 1794 for (i = 0; i < numIntFacets; i++) { 1795 const PetscInt cell0 = facetCells[2*(intFacetOffset + i) + 0]; 1796 const PetscInt cell1 = facetCells[2*(intFacetOffset + i) + 1]; 1797 PetscInt celli, cellj; 1798 1799 for (celli = 0; celli < patch->totalDofsPerCell; celli++) { 1800 const PetscInt row = dofsArray[(offset + cell0)*patch->totalDofsPerCell + celli]; 1801 if (row < 0) continue; 1802 for (cellj = 0; cellj < patch->totalDofsPerCell; cellj++) { 1803 const PetscInt col = dofsArray[(offset + cell1)*patch->totalDofsPerCell + cellj]; 1804 const PetscInt key = row*rsize + col; 1805 if (col < 0) continue; 1806 if (!PetscBTLookupSet(bt, key)) ++dnnz[row]; 1807 } 1808 } 1809 1810 for (celli = 0; celli < patch->totalDofsPerCell; celli++) { 1811 const PetscInt row = dofsArray[(offset + cell1)*patch->totalDofsPerCell + celli]; 1812 if (row < 0) continue; 1813 for (cellj = 0; cellj < patch->totalDofsPerCell; cellj++) { 1814 const PetscInt col = dofsArray[(offset + cell0)*patch->totalDofsPerCell + cellj]; 1815 const PetscInt key = row*rsize + col; 1816 if (col < 0) continue; 1817 if (!PetscBTLookupSet(bt, key)) ++dnnz[row]; 1818 } 1819 } 1820 } 1821 } 1822 ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 1823 ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr); 1824 ierr = PetscFree(dnnz);CHKERRQ(ierr); 1825 1826 ierr = PetscCalloc1(patch->totalDofsPerCell*patch->totalDofsPerCell, &zeroes);CHKERRQ(ierr); 1827 for (c = 0; c < ncell; ++c) { 1828 const PetscInt *idx = &dofsArray[(offset + c)*patch->totalDofsPerCell]; 1829 ierr = MatSetValues(*mat, patch->totalDofsPerCell, idx, patch->totalDofsPerCell, idx, zeroes, INSERT_VALUES);CHKERRQ(ierr); 1830 } 1831 ierr = MatGetLocalSize(*mat, &rows, NULL);CHKERRQ(ierr); 1832 for (i = 0; i < rows; ++i) { 1833 ierr = MatSetValues(*mat, 1, &i, 1, &i, zeroes, INSERT_VALUES);CHKERRQ(ierr); 1834 } 1835 1836 if (patch->usercomputeopintfacet) { 1837 const PetscInt *intFacetsArray = NULL; 1838 PetscInt i, numIntFacets, intFacetOffset; 1839 const PetscInt *facetCells = NULL; 1840 1841 ierr = PetscSectionGetDof(patch->intFacetCounts, point, &numIntFacets);CHKERRQ(ierr); 1842 ierr = PetscSectionGetOffset(patch->intFacetCounts, point, &intFacetOffset);CHKERRQ(ierr); 1843 ierr = ISGetIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 1844 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 1845 for (i = 0; i < numIntFacets; i++) { 1846 const PetscInt cell0 = facetCells[2*(intFacetOffset + i) + 0]; 1847 const PetscInt cell1 = facetCells[2*(intFacetOffset + i) + 1]; 1848 const PetscInt *cell0idx = &dofsArray[(offset + cell0)*patch->totalDofsPerCell]; 1849 const PetscInt *cell1idx = &dofsArray[(offset + cell1)*patch->totalDofsPerCell]; 1850 ierr = MatSetValues(*mat, patch->totalDofsPerCell, cell0idx, patch->totalDofsPerCell, cell1idx, zeroes, INSERT_VALUES);CHKERRQ(ierr); 1851 ierr = MatSetValues(*mat, patch->totalDofsPerCell, cell1idx, patch->totalDofsPerCell, cell0idx, zeroes, INSERT_VALUES);CHKERRQ(ierr); 1852 } 1853 } 1854 1855 ierr = MatAssemblyBegin(*mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1856 ierr = MatAssemblyEnd(*mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1857 1858 ierr = PetscFree(zeroes);CHKERRQ(ierr); 1859 1860 } else { /* rsize too big, use MATPREALLOCATOR */ 1861 Mat preallocator; 1862 PetscScalar* vals; 1863 1864 ierr = PetscCalloc1(patch->totalDofsPerCell*patch->totalDofsPerCell, &vals);CHKERRQ(ierr); 1865 ierr = MatCreate(PETSC_COMM_SELF, &preallocator);CHKERRQ(ierr); 1866 ierr = MatSetType(preallocator, MATPREALLOCATOR);CHKERRQ(ierr); 1867 ierr = MatSetSizes(preallocator, rsize, rsize, rsize, rsize);CHKERRQ(ierr); 1868 ierr = MatSetUp(preallocator);CHKERRQ(ierr); 1869 1870 for (c = 0; c < ncell; ++c) { 1871 const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell; 1872 ierr = MatSetValues(preallocator, patch->totalDofsPerCell, idx, patch->totalDofsPerCell, idx, vals, INSERT_VALUES);CHKERRQ(ierr); 1873 } 1874 1875 if (patch->usercomputeopintfacet) { 1876 const PetscInt *intFacetsArray = NULL; 1877 PetscInt i, numIntFacets, intFacetOffset; 1878 const PetscInt *facetCells = NULL; 1879 1880 ierr = PetscSectionGetDof(patch->intFacetCounts, point, &numIntFacets);CHKERRQ(ierr); 1881 ierr = PetscSectionGetOffset(patch->intFacetCounts, point, &intFacetOffset);CHKERRQ(ierr); 1882 ierr = ISGetIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 1883 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 1884 for (i = 0; i < numIntFacets; i++) { 1885 const PetscInt cell0 = facetCells[2*(intFacetOffset + i) + 0]; 1886 const PetscInt cell1 = facetCells[2*(intFacetOffset + i) + 1]; 1887 const PetscInt *cell0idx = &dofsArray[(offset + cell0)*patch->totalDofsPerCell]; 1888 const PetscInt *cell1idx = &dofsArray[(offset + cell1)*patch->totalDofsPerCell]; 1889 ierr = MatSetValues(preallocator, patch->totalDofsPerCell, cell0idx, patch->totalDofsPerCell, cell1idx, vals, INSERT_VALUES);CHKERRQ(ierr); 1890 ierr = MatSetValues(preallocator, patch->totalDofsPerCell, cell1idx, patch->totalDofsPerCell, cell0idx, vals, INSERT_VALUES);CHKERRQ(ierr); 1891 } 1892 } 1893 1894 ierr = PetscFree(vals);CHKERRQ(ierr); 1895 ierr = MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1896 ierr = MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1897 ierr = MatPreallocatorPreallocate(preallocator, PETSC_TRUE, *mat);CHKERRQ(ierr); 1898 ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 1899 } 1900 ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1901 if(withArtificial) { 1902 ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 1903 } else { 1904 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1905 } 1906 } 1907 ierr = MatSetUp(*mat);CHKERRQ(ierr); 1908 PetscFunctionReturn(0); 1909 } 1910 1911 static PetscErrorCode PCPatchComputeFunction_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Vec F, IS cellIS, PetscInt n, const PetscInt *l2p, const PetscInt *l2pWithAll, void *ctx) 1912 { 1913 PC_PATCH *patch = (PC_PATCH *) pc->data; 1914 DM dm; 1915 PetscSection s; 1916 const PetscInt *parray, *oarray; 1917 PetscInt Nf = patch->nsubspaces, Np, poff, p, f; 1918 PetscErrorCode ierr; 1919 1920 PetscFunctionBegin; 1921 if (patch->precomputeElementTensors) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Precomputing element tensors not implemented with DMPlex compute operator"); 1922 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 1923 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 1924 /* Set offset into patch */ 1925 ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr); 1926 ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr); 1927 ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr); 1928 ierr = ISGetIndices(patch->offs, &oarray);CHKERRQ(ierr); 1929 for (f = 0; f < Nf; ++f) { 1930 for (p = 0; p < Np; ++p) { 1931 const PetscInt point = parray[poff+p]; 1932 PetscInt dof; 1933 1934 ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr); 1935 ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr); 1936 if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);} 1937 else {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);} 1938 } 1939 } 1940 ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr); 1941 ierr = ISRestoreIndices(patch->offs, &oarray);CHKERRQ(ierr); 1942 if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);} 1943 ierr = DMPlexComputeResidual_Patch_Internal(pc->dm, patch->patchSection, cellIS, 0.0, x, NULL, F, ctx);CHKERRQ(ierr); 1944 PetscFunctionReturn(0); 1945 } 1946 1947 PetscErrorCode PCPatchComputeFunction_Internal(PC pc, Vec x, Vec F, PetscInt point) 1948 { 1949 PC_PATCH *patch = (PC_PATCH *) pc->data; 1950 const PetscInt *dofsArray; 1951 const PetscInt *dofsArrayWithAll; 1952 const PetscInt *cellsArray; 1953 PetscInt ncell, offset, pStart, pEnd; 1954 PetscErrorCode ierr; 1955 1956 PetscFunctionBegin; 1957 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1958 if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n"); 1959 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1960 ierr = ISGetIndices(patch->dofsWithAll, &dofsArrayWithAll);CHKERRQ(ierr); 1961 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1962 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1963 1964 point += pStart; 1965 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1966 1967 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1968 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1969 if (ncell <= 0) { 1970 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1971 PetscFunctionReturn(0); 1972 } 1973 PetscStackPush("PCPatch user callback"); 1974 /* Cannot reuse the same IS because the geometry info is being cached in it */ 1975 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr); 1976 ierr = patch->usercomputef(pc, point, x, F, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, 1977 dofsArrayWithAll + offset*patch->totalDofsPerCell, 1978 patch->usercomputefctx);CHKERRQ(ierr); 1979 PetscStackPop; 1980 ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr); 1981 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1982 ierr = ISRestoreIndices(patch->dofsWithAll, &dofsArrayWithAll);CHKERRQ(ierr); 1983 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1984 if (patch->viewMatrix) { 1985 char name[PETSC_MAX_PATH_LEN]; 1986 1987 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch vector for Point %D", point);CHKERRQ(ierr); 1988 ierr = PetscObjectSetName((PetscObject) F, name);CHKERRQ(ierr); 1989 ierr = ObjectView((PetscObject) F, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr); 1990 } 1991 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1992 PetscFunctionReturn(0); 1993 } 1994 1995 static PetscErrorCode PCPatchComputeOperator_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Mat J, IS cellIS, PetscInt n, const PetscInt *l2p, const PetscInt *l2pWithAll, void *ctx) 1996 { 1997 PC_PATCH *patch = (PC_PATCH *) pc->data; 1998 DM dm; 1999 PetscSection s; 2000 const PetscInt *parray, *oarray; 2001 PetscInt Nf = patch->nsubspaces, Np, poff, p, f; 2002 PetscErrorCode ierr; 2003 2004 PetscFunctionBegin; 2005 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 2006 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 2007 /* Set offset into patch */ 2008 ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr); 2009 ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr); 2010 ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr); 2011 ierr = ISGetIndices(patch->offs, &oarray);CHKERRQ(ierr); 2012 for (f = 0; f < Nf; ++f) { 2013 for (p = 0; p < Np; ++p) { 2014 const PetscInt point = parray[poff+p]; 2015 PetscInt dof; 2016 2017 ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr); 2018 ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr); 2019 if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);} 2020 else {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);} 2021 } 2022 } 2023 ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr); 2024 ierr = ISRestoreIndices(patch->offs, &oarray);CHKERRQ(ierr); 2025 if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);} 2026 /* TODO Shut off MatViewFromOptions() in MatAssemblyEnd() here */ 2027 ierr = DMPlexComputeJacobian_Patch_Internal(pc->dm, patch->patchSection, patch->patchSection, cellIS, 0.0, 0.0, x, NULL, J, J, ctx);CHKERRQ(ierr); 2028 PetscFunctionReturn(0); 2029 } 2030 2031 PetscErrorCode PCPatchComputeOperator_Internal(PC pc, Vec x, Mat mat, PetscInt point, PetscBool withArtificial) 2032 { 2033 PC_PATCH *patch = (PC_PATCH *) pc->data; 2034 const PetscInt *dofsArray; 2035 const PetscInt *dofsArrayWithAll = NULL; 2036 const PetscInt *cellsArray; 2037 PetscInt ncell, offset, pStart, pEnd, numIntFacets, intFacetOffset; 2038 PetscBool isNonlinear; 2039 PetscErrorCode ierr; 2040 2041 PetscFunctionBegin; 2042 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2043 isNonlinear = patch->isNonlinear; 2044 if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n"); 2045 if(withArtificial) { 2046 ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 2047 } else { 2048 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 2049 } 2050 if (isNonlinear) { 2051 ierr = ISGetIndices(patch->dofsWithAll, &dofsArrayWithAll);CHKERRQ(ierr); 2052 } 2053 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 2054 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 2055 2056 point += pStart; 2057 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 2058 2059 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 2060 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 2061 if (ncell <= 0) { 2062 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2063 PetscFunctionReturn(0); 2064 } 2065 if (patch->precomputeElementTensors) { 2066 PetscInt i; 2067 PetscInt ndof = patch->totalDofsPerCell; 2068 const PetscScalar *elementTensors; 2069 2070 ierr = VecGetArrayRead(patch->cellMats, &elementTensors);CHKERRQ(ierr); 2071 for (i = 0; i < ncell; i++) { 2072 const PetscInt cell = cellsArray[i + offset]; 2073 const PetscInt *idx = dofsArray + (offset + i)*ndof; 2074 const PetscScalar *v = elementTensors + patch->precomputedTensorLocations[cell]*ndof*ndof; 2075 ierr = MatSetValues(mat, ndof, idx, ndof, idx, v, ADD_VALUES);CHKERRQ(ierr); 2076 } 2077 ierr = VecRestoreArrayRead(patch->cellMats, &elementTensors);CHKERRQ(ierr); 2078 ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2079 ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2080 } else { 2081 PetscStackPush("PCPatch user callback"); 2082 /* Cannot reuse the same IS because the geometry info is being cached in it */ 2083 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr); 2084 ierr = patch->usercomputeop(pc, point, x, mat, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, dofsArrayWithAll ? dofsArrayWithAll + offset*patch->totalDofsPerCell : NULL, patch->usercomputeopctx);CHKERRQ(ierr); 2085 } 2086 if (patch->usercomputeopintfacet) { 2087 ierr = PetscSectionGetDof(patch->intFacetCounts, point, &numIntFacets);CHKERRQ(ierr); 2088 ierr = PetscSectionGetOffset(patch->intFacetCounts, point, &intFacetOffset);CHKERRQ(ierr); 2089 if (numIntFacets > 0) { 2090 /* For each interior facet, grab the two cells (in local numbering, and concatenate dof numberings for those cells) */ 2091 PetscInt *facetDofs = NULL, *facetDofsWithAll = NULL; 2092 const PetscInt *intFacetsArray = NULL; 2093 PetscInt idx = 0; 2094 PetscInt i, c, d; 2095 PetscInt fStart; 2096 DM dm; 2097 IS facetIS = NULL; 2098 const PetscInt *facetCells = NULL; 2099 2100 ierr = ISGetIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 2101 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 2102 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 2103 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, NULL);CHKERRQ(ierr); 2104 /* FIXME: Pull this malloc out. */ 2105 ierr = PetscMalloc1(2 * patch->totalDofsPerCell * numIntFacets, &facetDofs);CHKERRQ(ierr); 2106 if (dofsArrayWithAll) { 2107 ierr = PetscMalloc1(2 * patch->totalDofsPerCell * numIntFacets, &facetDofsWithAll);CHKERRQ(ierr); 2108 } 2109 if (patch->precomputeElementTensors) { 2110 PetscInt nFacetDof = 2*patch->totalDofsPerCell; 2111 const PetscScalar *elementTensors; 2112 2113 ierr = VecGetArrayRead(patch->intFacetMats, &elementTensors);CHKERRQ(ierr); 2114 2115 for (i = 0; i < numIntFacets; i++) { 2116 const PetscInt facet = intFacetsArray[i + intFacetOffset]; 2117 const PetscScalar *v = elementTensors + patch->precomputedIntFacetTensorLocations[facet - fStart]*nFacetDof*nFacetDof; 2118 idx = 0; 2119 /* 2120 * 0--1 2121 * |\-| 2122 * |+\| 2123 * 2--3 2124 * [0, 2, 3, 0, 1, 3] 2125 */ 2126 for (c = 0; c < 2; c++) { 2127 const PetscInt cell = facetCells[2*(intFacetOffset + i) + c]; 2128 for (d = 0; d < patch->totalDofsPerCell; d++) { 2129 facetDofs[idx] = dofsArray[(offset + cell)*patch->totalDofsPerCell + d]; 2130 idx++; 2131 } 2132 } 2133 ierr = MatSetValues(mat, nFacetDof, facetDofs, nFacetDof, facetDofs, v, ADD_VALUES);CHKERRQ(ierr); 2134 } 2135 ierr = VecRestoreArrayRead(patch->intFacetMats, &elementTensors);CHKERRQ(ierr); 2136 } else { 2137 /* 2138 * 0--1 2139 * |\-| 2140 * |+\| 2141 * 2--3 2142 * [0, 2, 3, 0, 1, 3] 2143 */ 2144 for (i = 0; i < numIntFacets; i++) { 2145 for (c = 0; c < 2; c++) { 2146 const PetscInt cell = facetCells[2*(intFacetOffset + i) + c]; 2147 for (d = 0; d < patch->totalDofsPerCell; d++) { 2148 facetDofs[idx] = dofsArray[(offset + cell)*patch->totalDofsPerCell + d]; 2149 if (dofsArrayWithAll) { 2150 facetDofsWithAll[idx] = dofsArrayWithAll[(offset + cell)*patch->totalDofsPerCell + d]; 2151 } 2152 idx++; 2153 } 2154 } 2155 } 2156 ierr = ISCreateGeneral(PETSC_COMM_SELF, numIntFacets, intFacetsArray + intFacetOffset, PETSC_USE_POINTER, &facetIS);CHKERRQ(ierr); 2157 ierr = patch->usercomputeopintfacet(pc, point, x, mat, facetIS, 2*numIntFacets*patch->totalDofsPerCell, facetDofs, facetDofsWithAll, patch->usercomputeopintfacetctx);CHKERRQ(ierr); 2158 ierr = ISDestroy(&facetIS);CHKERRQ(ierr); 2159 } 2160 ierr = ISRestoreIndices(patch->intFacetsToPatchCell, &facetCells);CHKERRQ(ierr); 2161 ierr = ISRestoreIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 2162 ierr = PetscFree(facetDofs);CHKERRQ(ierr); 2163 ierr = PetscFree(facetDofsWithAll);CHKERRQ(ierr); 2164 } 2165 } 2166 2167 ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2168 ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2169 2170 PetscStackPop; 2171 ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr); 2172 if(withArtificial) { 2173 ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 2174 } else { 2175 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 2176 } 2177 if (isNonlinear) { 2178 ierr = ISRestoreIndices(patch->dofsWithAll, &dofsArrayWithAll);CHKERRQ(ierr); 2179 } 2180 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 2181 if (patch->viewMatrix) { 2182 char name[PETSC_MAX_PATH_LEN]; 2183 2184 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch matrix for Point %D", point);CHKERRQ(ierr); 2185 ierr = PetscObjectSetName((PetscObject) mat, name);CHKERRQ(ierr); 2186 ierr = ObjectView((PetscObject) mat, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr); 2187 } 2188 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2189 PetscFunctionReturn(0); 2190 } 2191 2192 static PetscErrorCode MatSetValues_PCPatch_Private(Mat mat, PetscInt m, const PetscInt idxm[], 2193 PetscInt n, const PetscInt idxn[], const PetscScalar *v, InsertMode addv) 2194 { 2195 Vec data; 2196 PetscScalar *array; 2197 PetscInt bs, nz, i, j, cell; 2198 PetscErrorCode ierr; 2199 2200 ierr = MatShellGetContext(mat, &data);CHKERRQ(ierr); 2201 ierr = VecGetBlockSize(data, &bs);CHKERRQ(ierr); 2202 ierr = VecGetSize(data, &nz);CHKERRQ(ierr); 2203 ierr = VecGetArray(data, &array);CHKERRQ(ierr); 2204 if (m != n) SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Only for square insertion"); 2205 cell = (PetscInt)(idxm[0]/bs); // use the fact that this is called once per cell 2206 for (i = 0; i < m; i++) { 2207 //const PetscInt row = idxm[i]; 2208 if (idxm[i] != idxn[i]) SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Row and column indices must match!"); 2209 for (j = 0; j < n; j++) { 2210 const PetscScalar v_ = v[i*bs + j]; 2211 /* Indexing is special to the data structure we have! */ 2212 if (addv == INSERT_VALUES) { 2213 array[cell*bs*bs + i*bs + j] = v_; 2214 } else { 2215 array[cell*bs*bs + i*bs + j] += v_; 2216 } 2217 } 2218 } 2219 ierr = VecRestoreArray(data, &array);CHKERRQ(ierr); 2220 PetscFunctionReturn(0); 2221 } 2222 2223 static PetscErrorCode PCPatchPrecomputePatchTensors_Private(PC pc) 2224 { 2225 PC_PATCH *patch = (PC_PATCH *)pc->data; 2226 const PetscInt *cellsArray; 2227 PetscInt ncell, offset; 2228 const PetscInt *dofMapArray; 2229 PetscInt i, j; 2230 IS dofMap; 2231 IS cellIS; 2232 const PetscInt ndof = patch->totalDofsPerCell; 2233 PetscErrorCode ierr; 2234 Mat vecMat; 2235 PetscInt cStart, cEnd; 2236 DM dm, plex; 2237 2238 2239 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2240 2241 if (!patch->allCells) { 2242 PetscHSetI cells; 2243 PetscHashIter hi; 2244 PetscInt pStart, pEnd; 2245 PetscInt *allCells = NULL; 2246 ierr = PetscHSetICreate(&cells);CHKERRQ(ierr); 2247 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 2248 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 2249 for (i = pStart; i < pEnd; i++) { 2250 ierr = PetscSectionGetDof(patch->cellCounts, i, &ncell);CHKERRQ(ierr); 2251 ierr = PetscSectionGetOffset(patch->cellCounts, i, &offset);CHKERRQ(ierr); 2252 if (ncell <= 0) continue; 2253 for (j = 0; j < ncell; j++) { 2254 PetscHSetIAdd(cells, cellsArray[offset + j]);CHKERRQ(ierr); 2255 } 2256 } 2257 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 2258 ierr = PetscHSetIGetSize(cells, &ncell);CHKERRQ(ierr); 2259 ierr = PetscMalloc1(ncell, &allCells);CHKERRQ(ierr); 2260 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 2261 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 2262 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2263 ierr = PetscMalloc1(cEnd-cStart, &patch->precomputedTensorLocations);CHKERRQ(ierr); 2264 i = 0; 2265 PetscHashIterBegin(cells, hi); 2266 while (!PetscHashIterAtEnd(cells, hi)) { 2267 PetscHashIterGetKey(cells, hi, allCells[i]); 2268 patch->precomputedTensorLocations[allCells[i]] = i; 2269 PetscHashIterNext(cells, hi); 2270 i++; 2271 } 2272 ierr = PetscHSetIDestroy(&cells);CHKERRQ(ierr); 2273 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, allCells, PETSC_OWN_POINTER, &patch->allCells);CHKERRQ(ierr); 2274 } 2275 ierr = ISGetSize(patch->allCells, &ncell);CHKERRQ(ierr); 2276 if (!patch->cellMats) { 2277 ierr = VecCreateSeq(PETSC_COMM_SELF, ncell*ndof*ndof, &patch->cellMats);CHKERRQ(ierr); 2278 ierr = VecSetBlockSize(patch->cellMats, ndof);CHKERRQ(ierr); 2279 } 2280 ierr = VecSet(patch->cellMats, 0);CHKERRQ(ierr); 2281 2282 ierr = MatCreateShell(PETSC_COMM_SELF, ncell*ndof, ncell*ndof, ncell*ndof, ncell*ndof, 2283 (void*)patch->cellMats, &vecMat); 2284 ierr = MatShellSetOperation(vecMat, MATOP_SET_VALUES, (void(*)(void))&MatSetValues_PCPatch_Private);CHKERRQ(ierr); 2285 ierr = ISGetSize(patch->allCells, &ncell);CHKERRQ(ierr); 2286 ierr = ISCreateStride(PETSC_COMM_SELF, ndof*ncell, 0, 1, &dofMap);CHKERRQ(ierr); 2287 ierr = ISGetIndices(dofMap, &dofMapArray);CHKERRQ(ierr); 2288 ierr = ISGetIndices(patch->allCells, &cellsArray);CHKERRQ(ierr); 2289 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray, PETSC_USE_POINTER, &cellIS);CHKERRQ(ierr); 2290 PetscStackPush("PCPatch user callback"); 2291 /* TODO: Fix for DMPlex compute op, this bypasses a lot of the machinery and just assembles every element tensor. */ 2292 ierr = patch->usercomputeop(pc, -1, NULL, vecMat, cellIS, ndof*ncell, dofMapArray, NULL, patch->usercomputeopctx);CHKERRQ(ierr); 2293 PetscStackPop; 2294 ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 2295 ierr = MatDestroy(&vecMat);CHKERRQ(ierr); 2296 ierr = ISRestoreIndices(patch->allCells, &cellsArray);CHKERRQ(ierr); 2297 ierr = ISRestoreIndices(dofMap, &dofMapArray);CHKERRQ(ierr); 2298 ierr = ISDestroy(&dofMap);CHKERRQ(ierr); 2299 2300 if (patch->usercomputeopintfacet) { 2301 PetscInt nIntFacets; 2302 IS intFacetsIS; 2303 const PetscInt *intFacetsArray = NULL; 2304 if (!patch->allIntFacets) { 2305 PetscHSetI facets; 2306 PetscHashIter hi; 2307 PetscInt pStart, pEnd, fStart, fEnd; 2308 PetscInt *allIntFacets = NULL; 2309 ierr = PetscHSetICreate(&facets);CHKERRQ(ierr); 2310 ierr = ISGetIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 2311 ierr = PetscSectionGetChart(patch->intFacetCounts, &pStart, &pEnd);CHKERRQ(ierr); 2312 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 2313 for (i = pStart; i < pEnd; i++) { 2314 ierr = PetscSectionGetDof(patch->intFacetCounts, i, &nIntFacets);CHKERRQ(ierr); 2315 ierr = PetscSectionGetOffset(patch->intFacetCounts, i, &offset);CHKERRQ(ierr); 2316 if (nIntFacets <= 0) continue; 2317 for (j = 0; j < nIntFacets; j++) { 2318 PetscHSetIAdd(facets, intFacetsArray[offset + j]);CHKERRQ(ierr); 2319 } 2320 } 2321 ierr = ISRestoreIndices(patch->intFacets, &intFacetsArray);CHKERRQ(ierr); 2322 ierr = PetscHSetIGetSize(facets, &nIntFacets);CHKERRQ(ierr); 2323 ierr = PetscMalloc1(nIntFacets, &allIntFacets);CHKERRQ(ierr); 2324 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 2325 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 2326 ierr = PetscMalloc1(fEnd-fStart, &patch->precomputedIntFacetTensorLocations);CHKERRQ(ierr); 2327 i = 0; 2328 PetscHashIterBegin(facets, hi); 2329 while (!PetscHashIterAtEnd(facets, hi)) { 2330 PetscHashIterGetKey(facets, hi, allIntFacets[i]); 2331 patch->precomputedIntFacetTensorLocations[allIntFacets[i] - fStart] = i; 2332 PetscHashIterNext(facets, hi); 2333 i++; 2334 } 2335 ierr = PetscHSetIDestroy(&facets);CHKERRQ(ierr); 2336 ierr = ISCreateGeneral(PETSC_COMM_SELF, nIntFacets, allIntFacets, PETSC_OWN_POINTER, &patch->allIntFacets);CHKERRQ(ierr); 2337 } 2338 ierr = ISGetSize(patch->allIntFacets, &nIntFacets);CHKERRQ(ierr); 2339 if (!patch->intFacetMats) { 2340 ierr = VecCreateSeq(PETSC_COMM_SELF, nIntFacets*ndof*ndof*4, &patch->intFacetMats);CHKERRQ(ierr); 2341 ierr = VecSetBlockSize(patch->intFacetMats, ndof*2);CHKERRQ(ierr); 2342 } 2343 ierr = VecSet(patch->intFacetMats, 0);CHKERRQ(ierr); 2344 2345 ierr = MatCreateShell(PETSC_COMM_SELF, nIntFacets*ndof*2, nIntFacets*ndof*2, nIntFacets*ndof*2, nIntFacets*ndof*2, 2346 (void*)patch->intFacetMats, &vecMat); 2347 ierr = MatShellSetOperation(vecMat, MATOP_SET_VALUES, (void(*)(void))&MatSetValues_PCPatch_Private);CHKERRQ(ierr); 2348 ierr = ISCreateStride(PETSC_COMM_SELF, 2*ndof*nIntFacets, 0, 1, &dofMap);CHKERRQ(ierr); 2349 ierr = ISGetIndices(dofMap, &dofMapArray);CHKERRQ(ierr); 2350 ierr = ISGetIndices(patch->allIntFacets, &intFacetsArray);CHKERRQ(ierr); 2351 ierr = ISCreateGeneral(PETSC_COMM_SELF, nIntFacets, intFacetsArray, PETSC_USE_POINTER, &intFacetsIS);CHKERRQ(ierr); 2352 PetscStackPush("PCPatch user callback (interior facets)"); 2353 /* TODO: Fix for DMPlex compute op, this bypasses a lot of the machinery and just assembles every element tensor. */ 2354 ierr = patch->usercomputeopintfacet(pc, -1, NULL, vecMat, intFacetsIS, 2*ndof*nIntFacets, dofMapArray, NULL, patch->usercomputeopintfacetctx);CHKERRQ(ierr); 2355 PetscStackPop; 2356 ierr = ISDestroy(&intFacetsIS);CHKERRQ(ierr); 2357 ierr = MatDestroy(&vecMat);CHKERRQ(ierr); 2358 ierr = ISRestoreIndices(patch->allIntFacets, &intFacetsArray);CHKERRQ(ierr); 2359 ierr = ISRestoreIndices(dofMap, &dofMapArray);CHKERRQ(ierr); 2360 ierr = ISDestroy(&dofMap);CHKERRQ(ierr); 2361 } 2362 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 2363 2364 PetscFunctionReturn(0); 2365 } 2366 2367 PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat, PatchScatterType scattertype) 2368 { 2369 PC_PATCH *patch = (PC_PATCH *) pc->data; 2370 const PetscScalar *xArray = NULL; 2371 PetscScalar *yArray = NULL; 2372 const PetscInt *gtolArray = NULL; 2373 PetscInt dof, offset, lidx; 2374 PetscErrorCode ierr; 2375 2376 PetscFunctionBeginHot; 2377 ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 2378 ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr); 2379 ierr = VecGetArray(y, &yArray);CHKERRQ(ierr); 2380 if (scattertype == SCATTER_WITHARTIFICIAL) { 2381 ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr); 2382 ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offset);CHKERRQ(ierr); 2383 ierr = ISGetIndices(patch->gtolWithArtificial, >olArray);CHKERRQ(ierr); 2384 } else if (scattertype == SCATTER_WITHALL) { 2385 ierr = PetscSectionGetDof(patch->gtolCountsWithAll, p, &dof);CHKERRQ(ierr); 2386 ierr = PetscSectionGetOffset(patch->gtolCountsWithAll, p, &offset);CHKERRQ(ierr); 2387 ierr = ISGetIndices(patch->gtolWithAll, >olArray);CHKERRQ(ierr); 2388 } else { 2389 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 2390 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 2391 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2392 } 2393 if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n"); 2394 if (mode == ADD_VALUES && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n"); 2395 for (lidx = 0; lidx < dof; ++lidx) { 2396 const PetscInt gidx = gtolArray[offset+lidx]; 2397 2398 if (mode == INSERT_VALUES) yArray[lidx] = xArray[gidx]; /* Forward */ 2399 else yArray[gidx] += xArray[lidx]; /* Reverse */ 2400 } 2401 if (scattertype == SCATTER_WITHARTIFICIAL) { 2402 ierr = ISRestoreIndices(patch->gtolWithArtificial, >olArray);CHKERRQ(ierr); 2403 } else if (scattertype == SCATTER_WITHALL) { 2404 ierr = ISRestoreIndices(patch->gtolWithAll, >olArray);CHKERRQ(ierr); 2405 } else { 2406 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2407 } 2408 ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr); 2409 ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr); 2410 ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 2411 PetscFunctionReturn(0); 2412 } 2413 2414 static PetscErrorCode PCSetUp_PATCH_Linear(PC pc) 2415 { 2416 PC_PATCH *patch = (PC_PATCH *) pc->data; 2417 const char *prefix; 2418 PetscInt i; 2419 PetscErrorCode ierr; 2420 2421 PetscFunctionBegin; 2422 if (!pc->setupcalled) { 2423 ierr = PetscMalloc1(patch->npatch, &patch->solver);CHKERRQ(ierr); 2424 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 2425 for (i = 0; i < patch->npatch; ++i) { 2426 KSP ksp; 2427 PC subpc; 2428 2429 ierr = KSPCreate(PETSC_COMM_SELF, &ksp);CHKERRQ(ierr); 2430 ierr = KSPSetErrorIfNotConverged(ksp, pc->erroriffailure);CHKERRQ(ierr); 2431 ierr = KSPSetOptionsPrefix(ksp, prefix);CHKERRQ(ierr); 2432 ierr = KSPAppendOptionsPrefix(ksp, "sub_");CHKERRQ(ierr); 2433 ierr = PetscObjectIncrementTabLevel((PetscObject) ksp, (PetscObject) pc, 1);CHKERRQ(ierr); 2434 ierr = KSPGetPC(ksp, &subpc);CHKERRQ(ierr); 2435 ierr = PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1);CHKERRQ(ierr); 2436 ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) ksp);CHKERRQ(ierr); 2437 patch->solver[i] = (PetscObject) ksp; 2438 } 2439 } 2440 if (patch->save_operators) { 2441 if (patch->precomputeElementTensors) { 2442 ierr = PCPatchPrecomputePatchTensors_Private(pc);CHKERRQ(ierr); 2443 } 2444 for (i = 0; i < patch->npatch; ++i) { 2445 ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr); 2446 ierr = PCPatchComputeOperator_Internal(pc, NULL, patch->mat[i], i, PETSC_FALSE);CHKERRQ(ierr); 2447 ierr = KSPSetOperators((KSP) patch->solver[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr); 2448 } 2449 } 2450 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2451 for (i = 0; i < patch->npatch; ++i) { 2452 /* Instead of padding patch->patchUpdate with zeros to get */ 2453 /* patch->patchUpdateWithArtificial and then multiplying with the matrix, */ 2454 /* just get rid of the columns that correspond to the dofs with */ 2455 /* artificial bcs. That's of course fairly inefficient, hopefully we */ 2456 /* can just assemble the rectangular matrix in the first place. */ 2457 Mat matSquare; 2458 IS rowis; 2459 PetscInt dof; 2460 2461 ierr = MatGetSize(patch->mat[i], &dof, NULL);CHKERRQ(ierr); 2462 if (dof == 0) { 2463 patch->matWithArtificial[i] = NULL; 2464 continue; 2465 } 2466 2467 ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr); 2468 ierr = MatZeroEntries(matSquare);CHKERRQ(ierr); 2469 ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr); 2470 2471 ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr); 2472 ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr); 2473 if(pc->setupcalled) { 2474 ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_REUSE_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr); 2475 } else { 2476 ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr); 2477 } 2478 ierr = ISDestroy(&rowis); CHKERRQ(ierr); 2479 ierr = MatDestroy(&matSquare);CHKERRQ(ierr); 2480 } 2481 } 2482 PetscFunctionReturn(0); 2483 } 2484 2485 static PetscErrorCode PCSetUp_PATCH(PC pc) 2486 { 2487 PC_PATCH *patch = (PC_PATCH *) pc->data; 2488 PetscInt i; 2489 PetscBool isNonlinear; 2490 PetscErrorCode ierr; 2491 2492 PetscFunctionBegin; 2493 if (!pc->setupcalled) { 2494 PetscInt pStart, pEnd, p; 2495 PetscInt localSize; 2496 2497 ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 2498 2499 isNonlinear = patch->isNonlinear; 2500 if (!patch->nsubspaces) { 2501 DM dm; 2502 PetscSection s; 2503 PetscInt cStart, cEnd, c, Nf, f, numGlobalBcs = 0, *globalBcs, *Nb, totNb = 0, **cellDofs; 2504 2505 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 2506 if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Must set DM for PCPATCH or call PCPatchSetDiscretisationInfo()"); 2507 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 2508 ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr); 2509 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 2510 for (p = pStart; p < pEnd; ++p) { 2511 PetscInt cdof; 2512 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 2513 numGlobalBcs += cdof; 2514 } 2515 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 2516 ierr = PetscMalloc3(Nf, &Nb, Nf, &cellDofs, numGlobalBcs, &globalBcs);CHKERRQ(ierr); 2517 for (f = 0; f < Nf; ++f) { 2518 PetscFE fe; 2519 PetscDualSpace sp; 2520 PetscInt cdoff = 0; 2521 2522 ierr = DMGetField(dm, f, NULL, (PetscObject *) &fe);CHKERRQ(ierr); 2523 /* ierr = PetscFEGetNumComponents(fe, &Nc[f]);CHKERRQ(ierr); */ 2524 ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 2525 ierr = PetscDualSpaceGetDimension(sp, &Nb[f]);CHKERRQ(ierr); 2526 totNb += Nb[f]; 2527 2528 ierr = PetscMalloc1((cEnd-cStart)*Nb[f], &cellDofs[f]);CHKERRQ(ierr); 2529 for (c = cStart; c < cEnd; ++c) { 2530 PetscInt *closure = NULL; 2531 PetscInt clSize = 0, cl; 2532 2533 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 2534 for (cl = 0; cl < clSize*2; cl += 2) { 2535 const PetscInt p = closure[cl]; 2536 PetscInt fdof, d, foff; 2537 2538 ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 2539 ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr); 2540 for (d = 0; d < fdof; ++d, ++cdoff) cellDofs[f][cdoff] = foff + d; 2541 } 2542 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 2543 } 2544 if (cdoff != (cEnd-cStart)*Nb[f]) SETERRQ4(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_SIZ, "Total number of cellDofs %D for field %D should be Nc (%D) * cellDof (%D)", cdoff, f, cEnd-cStart, Nb[f]); 2545 } 2546 numGlobalBcs = 0; 2547 for (p = pStart; p < pEnd; ++p) { 2548 const PetscInt *ind; 2549 PetscInt off, cdof, d; 2550 2551 ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr); 2552 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 2553 ierr = PetscSectionGetConstraintIndices(s, p, &ind);CHKERRQ(ierr); 2554 for (d = 0; d < cdof; ++d) globalBcs[numGlobalBcs++] = off + ind[d]; 2555 } 2556 2557 ierr = PCPatchSetDiscretisationInfoCombined(pc, dm, Nb, (const PetscInt **) cellDofs, numGlobalBcs, globalBcs, numGlobalBcs, globalBcs);CHKERRQ(ierr); 2558 for (f = 0; f < Nf; ++f) { 2559 ierr = PetscFree(cellDofs[f]);CHKERRQ(ierr); 2560 } 2561 ierr = PetscFree3(Nb, cellDofs, globalBcs);CHKERRQ(ierr); 2562 ierr = PCPatchSetComputeFunction(pc, PCPatchComputeFunction_DMPlex_Private, NULL);CHKERRQ(ierr); 2563 ierr = PCPatchSetComputeOperator(pc, PCPatchComputeOperator_DMPlex_Private, NULL);CHKERRQ(ierr); 2564 } 2565 2566 localSize = patch->subspaceOffsets[patch->nsubspaces]; 2567 ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localRHS);CHKERRQ(ierr); 2568 ierr = VecSetUp(patch->localRHS);CHKERRQ(ierr); 2569 ierr = VecDuplicate(patch->localRHS, &patch->localUpdate);CHKERRQ(ierr); 2570 ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr); 2571 ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr); 2572 2573 /* OK, now build the work vectors */ 2574 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr); 2575 ierr = PetscMalloc1(patch->npatch, &patch->patchRHS);CHKERRQ(ierr); 2576 ierr = PetscMalloc1(patch->npatch, &patch->patchUpdate);CHKERRQ(ierr); 2577 2578 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2579 ierr = PetscMalloc1(patch->npatch, &patch->patchRHSWithArtificial);CHKERRQ(ierr); 2580 ierr = PetscMalloc1(patch->npatch, &patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr); 2581 } 2582 if (isNonlinear) { 2583 ierr = PetscMalloc1(patch->npatch, &patch->dofMappingWithoutToWithAll);CHKERRQ(ierr); 2584 } 2585 for (p = pStart; p < pEnd; ++p) { 2586 PetscInt dof; 2587 2588 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 2589 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHS[p-pStart]);CHKERRQ(ierr); 2590 ierr = VecSetUp(patch->patchRHS[p-pStart]);CHKERRQ(ierr); 2591 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchUpdate[p-pStart]);CHKERRQ(ierr); 2592 ierr = VecSetUp(patch->patchUpdate[p-pStart]);CHKERRQ(ierr); 2593 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2594 const PetscInt *gtolArray, *gtolArrayWithArtificial = NULL; 2595 PetscInt numPatchDofs, offset; 2596 PetscInt numPatchDofsWithArtificial, offsetWithArtificial; 2597 PetscInt dofWithoutArtificialCounter = 0; 2598 PetscInt *patchWithoutArtificialToWithArtificialArray; 2599 2600 ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr); 2601 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr); 2602 ierr = VecSetUp(patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr); 2603 2604 /* Now build the mapping that for a dof in a patch WITHOUT dofs that have artificial bcs gives the */ 2605 /* the index in the patch with all dofs */ 2606 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2607 2608 ierr = PetscSectionGetDof(patch->gtolCounts, p, &numPatchDofs);CHKERRQ(ierr); 2609 if (numPatchDofs == 0) { 2610 patch->dofMappingWithoutToWithArtificial[p-pStart] = NULL; 2611 continue; 2612 } 2613 2614 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 2615 ierr = ISGetIndices(patch->gtolWithArtificial, >olArrayWithArtificial);CHKERRQ(ierr); 2616 ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &numPatchDofsWithArtificial);CHKERRQ(ierr); 2617 ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offsetWithArtificial);CHKERRQ(ierr); 2618 2619 ierr = PetscMalloc1(numPatchDofs, &patchWithoutArtificialToWithArtificialArray);CHKERRQ(ierr); 2620 for (i=0; i<numPatchDofsWithArtificial; i++) { 2621 if (gtolArrayWithArtificial[i+offsetWithArtificial] == gtolArray[offset+dofWithoutArtificialCounter]) { 2622 patchWithoutArtificialToWithArtificialArray[dofWithoutArtificialCounter] = i; 2623 dofWithoutArtificialCounter++; 2624 if (dofWithoutArtificialCounter == numPatchDofs) 2625 break; 2626 } 2627 } 2628 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPatchDofs, patchWithoutArtificialToWithArtificialArray, PETSC_OWN_POINTER, &patch->dofMappingWithoutToWithArtificial[p-pStart]);CHKERRQ(ierr); 2629 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2630 ierr = ISRestoreIndices(patch->gtolWithArtificial, >olArrayWithArtificial);CHKERRQ(ierr); 2631 } 2632 if (isNonlinear) { 2633 const PetscInt *gtolArray, *gtolArrayWithAll = NULL; 2634 PetscInt numPatchDofs, offset; 2635 PetscInt numPatchDofsWithAll, offsetWithAll; 2636 PetscInt dofWithoutAllCounter = 0; 2637 PetscInt *patchWithoutAllToWithAllArray; 2638 2639 /* Now build the mapping that for a dof in a patch WITHOUT dofs that have artificial bcs gives the */ 2640 /* the index in the patch with all dofs */ 2641 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2642 2643 ierr = PetscSectionGetDof(patch->gtolCounts, p, &numPatchDofs);CHKERRQ(ierr); 2644 if (numPatchDofs == 0) { 2645 patch->dofMappingWithoutToWithAll[p-pStart] = NULL; 2646 continue; 2647 } 2648 2649 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 2650 ierr = ISGetIndices(patch->gtolWithAll, >olArrayWithAll);CHKERRQ(ierr); 2651 ierr = PetscSectionGetDof(patch->gtolCountsWithAll, p, &numPatchDofsWithAll);CHKERRQ(ierr); 2652 ierr = PetscSectionGetOffset(patch->gtolCountsWithAll, p, &offsetWithAll);CHKERRQ(ierr); 2653 2654 ierr = PetscMalloc1(numPatchDofs, &patchWithoutAllToWithAllArray);CHKERRQ(ierr); 2655 2656 for (i=0; i<numPatchDofsWithAll; i++) { 2657 if (gtolArrayWithAll[i+offsetWithAll] == gtolArray[offset+dofWithoutAllCounter]) { 2658 patchWithoutAllToWithAllArray[dofWithoutAllCounter] = i; 2659 dofWithoutAllCounter++; 2660 if (dofWithoutAllCounter == numPatchDofs) 2661 break; 2662 } 2663 } 2664 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPatchDofs, patchWithoutAllToWithAllArray, PETSC_OWN_POINTER, &patch->dofMappingWithoutToWithAll[p-pStart]);CHKERRQ(ierr); 2665 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 2666 ierr = ISRestoreIndices(patch->gtolWithAll, >olArrayWithAll);CHKERRQ(ierr); 2667 } 2668 } 2669 if (patch->save_operators) { 2670 ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr); 2671 for (i = 0; i < patch->npatch; ++i) { 2672 ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i], PETSC_FALSE);CHKERRQ(ierr); 2673 } 2674 } 2675 ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 2676 2677 /* If desired, calculate weights for dof multiplicity */ 2678 if (patch->partition_of_unity) { 2679 PetscScalar *input = NULL; 2680 PetscScalar *output = NULL; 2681 Vec global; 2682 2683 ierr = VecDuplicate(patch->localRHS, &patch->dof_weights);CHKERRQ(ierr); 2684 if(patch->local_composition_type == PC_COMPOSITE_ADDITIVE) { 2685 for (i = 0; i < patch->npatch; ++i) { 2686 PetscInt dof; 2687 2688 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr); 2689 if (dof <= 0) continue; 2690 ierr = VecSet(patch->patchRHS[i], 1.0);CHKERRQ(ierr); 2691 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchRHS[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE, SCATTER_INTERIOR);CHKERRQ(ierr); 2692 } 2693 } else { 2694 /* multiplicative is actually only locally multiplicative and globally additive. need the pou where the mesh decomposition overlaps */ 2695 ierr = VecSet(patch->dof_weights, 1.0);CHKERRQ(ierr); 2696 } 2697 2698 VecDuplicate(patch->dof_weights, &global); 2699 VecSet(global, 0.); 2700 2701 ierr = VecGetArray(patch->dof_weights, &input);CHKERRQ(ierr); 2702 ierr = VecGetArray(global, &output);CHKERRQ(ierr); 2703 ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr); 2704 ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr); 2705 ierr = VecRestoreArray(patch->dof_weights, &input);CHKERRQ(ierr); 2706 ierr = VecRestoreArray(global, &output);CHKERRQ(ierr); 2707 2708 ierr = VecReciprocal(global);CHKERRQ(ierr); 2709 2710 ierr = VecGetArray(patch->dof_weights, &output);CHKERRQ(ierr); 2711 ierr = VecGetArray(global, &input);CHKERRQ(ierr); 2712 ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr); 2713 ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr); 2714 ierr = VecRestoreArray(patch->dof_weights, &output);CHKERRQ(ierr); 2715 ierr = VecRestoreArray(global, &input);CHKERRQ(ierr); 2716 ierr = VecDestroy(&global);CHKERRQ(ierr); 2717 } 2718 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE && patch->save_operators) { 2719 ierr = PetscMalloc1(patch->npatch, &patch->matWithArtificial);CHKERRQ(ierr); 2720 } 2721 } 2722 ierr = (*patch->setupsolver)(pc);CHKERRQ(ierr); 2723 PetscFunctionReturn(0); 2724 } 2725 2726 static PetscErrorCode PCApply_PATCH_Linear(PC pc, PetscInt i, Vec x, Vec y) 2727 { 2728 PC_PATCH *patch = (PC_PATCH *) pc->data; 2729 KSP ksp = (KSP) patch->solver[i]; 2730 PetscErrorCode ierr; 2731 2732 PetscFunctionBegin; 2733 if (!patch->save_operators) { 2734 Mat mat; 2735 2736 ierr = PCPatchCreateMatrix_Private(pc, i, &mat, PETSC_FALSE);CHKERRQ(ierr); 2737 /* Populate operator here. */ 2738 ierr = PCPatchComputeOperator_Internal(pc, NULL, mat, i, PETSC_FALSE);CHKERRQ(ierr); 2739 ierr = KSPSetOperators(ksp, mat, mat);CHKERRQ(ierr); 2740 /* Drop reference so the KSPSetOperators below will blow it away. */ 2741 ierr = MatDestroy(&mat);CHKERRQ(ierr); 2742 } 2743 ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 2744 if (!ksp->setfromoptionscalled) { 2745 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 2746 } 2747 ierr = KSPSolve(ksp, x, y);CHKERRQ(ierr); 2748 ierr = KSPCheckSolve(ksp, pc, y);CHKERRQ(ierr); 2749 ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 2750 if (!patch->save_operators) { 2751 PC pc; 2752 ierr = KSPSetOperators(ksp, NULL, NULL);CHKERRQ(ierr); 2753 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 2754 /* Destroy PC context too, otherwise the factored matrix hangs around. */ 2755 ierr = PCReset(pc);CHKERRQ(ierr); 2756 } 2757 PetscFunctionReturn(0); 2758 } 2759 2760 static PetscErrorCode PCUpdateMultiplicative_PATCH_Linear(PC pc, PetscInt i, PetscInt pStart) 2761 { 2762 PC_PATCH *patch = (PC_PATCH *) pc->data; 2763 Mat multMat; 2764 PetscErrorCode ierr; 2765 2766 PetscFunctionBegin; 2767 2768 if (patch->save_operators) { 2769 multMat = patch->matWithArtificial[i]; 2770 } else { 2771 /*Very inefficient, hopefully we can just assemble the rectangular matrix in the first place.*/ 2772 Mat matSquare; 2773 PetscInt dof; 2774 IS rowis; 2775 ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr); 2776 ierr = MatZeroEntries(matSquare);CHKERRQ(ierr); 2777 ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr); 2778 ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr); 2779 ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr); 2780 ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &multMat); CHKERRQ(ierr); 2781 ierr = MatDestroy(&matSquare);CHKERRQ(ierr); 2782 ierr = ISDestroy(&rowis); CHKERRQ(ierr); 2783 } 2784 ierr = MatMult(multMat, patch->patchUpdate[i], patch->patchRHSWithArtificial[i]); CHKERRQ(ierr); 2785 ierr = VecScale(patch->patchRHSWithArtificial[i], -1.0); CHKERRQ(ierr); 2786 ierr = PCPatch_ScatterLocal_Private(pc, i + pStart, patch->patchRHSWithArtificial[i], patch->localRHS, ADD_VALUES, SCATTER_REVERSE, SCATTER_WITHARTIFICIAL); CHKERRQ(ierr); 2787 if (!patch->save_operators) { 2788 ierr = MatDestroy(&multMat); CHKERRQ(ierr); 2789 } 2790 PetscFunctionReturn(0); 2791 } 2792 2793 static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y) 2794 { 2795 PC_PATCH *patch = (PC_PATCH *) pc->data; 2796 const PetscScalar *globalRHS = NULL; 2797 PetscScalar *localRHS = NULL; 2798 PetscScalar *globalUpdate = NULL; 2799 const PetscInt *bcNodes = NULL; 2800 PetscInt nsweep = patch->symmetrise_sweep ? 2 : 1; 2801 PetscInt start[2] = {0, 0}; 2802 PetscInt end[2] = {-1, -1}; 2803 const PetscInt inc[2] = {1, -1}; 2804 const PetscScalar *localUpdate; 2805 const PetscInt *iterationSet; 2806 PetscInt pStart, numBcs, n, sweep, bc, j; 2807 PetscErrorCode ierr; 2808 2809 PetscFunctionBegin; 2810 ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 2811 ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr); 2812 /* start, end, inc have 2 entries to manage a second backward sweep if we symmetrize */ 2813 end[0] = patch->npatch; 2814 start[1] = patch->npatch-1; 2815 if (patch->user_patches) { 2816 ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr); 2817 start[1] = end[0] - 1; 2818 ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr); 2819 } 2820 /* Scatter from global space into overlapped local spaces */ 2821 ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr); 2822 ierr = VecGetArray(patch->localRHS, &localRHS);CHKERRQ(ierr); 2823 ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr); 2824 ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr); 2825 ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr); 2826 ierr = VecRestoreArray(patch->localRHS, &localRHS);CHKERRQ(ierr); 2827 2828 ierr = VecSet(patch->localUpdate, 0.0);CHKERRQ(ierr); 2829 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr); 2830 for (sweep = 0; sweep < nsweep; sweep++) { 2831 for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) { 2832 PetscInt i = patch->user_patches ? iterationSet[j] : j; 2833 PetscInt start, len; 2834 2835 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr); 2836 ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr); 2837 /* TODO: Squash out these guys in the setup as well. */ 2838 if (len <= 0) continue; 2839 /* TODO: Do we need different scatters for X and Y? */ 2840 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localRHS, patch->patchRHS[i], INSERT_VALUES, SCATTER_FORWARD, SCATTER_INTERIOR);CHKERRQ(ierr); 2841 ierr = (*patch->applysolver)(pc, i, patch->patchRHS[i], patch->patchUpdate[i]);CHKERRQ(ierr); 2842 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchUpdate[i], patch->localUpdate, ADD_VALUES, SCATTER_REVERSE, SCATTER_INTERIOR);CHKERRQ(ierr); 2843 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2844 ierr = (*patch->updatemultiplicative)(pc, i, pStart);CHKERRQ(ierr); 2845 } 2846 } 2847 } 2848 if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);} 2849 /* XXX: should we do this on the global vector? */ 2850 if (patch->partition_of_unity) { 2851 ierr = VecPointwiseMult(patch->localUpdate, patch->localUpdate, patch->dof_weights);CHKERRQ(ierr); 2852 } 2853 /* Now patch->localUpdate contains the solution of the patch solves, so we need to combine them all. */ 2854 ierr = VecSet(y, 0.0);CHKERRQ(ierr); 2855 ierr = VecGetArray(y, &globalUpdate);CHKERRQ(ierr); 2856 ierr = VecGetArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr); 2857 ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr); 2858 ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr); 2859 ierr = VecRestoreArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr); 2860 2861 /* Now we need to send the global BC values through */ 2862 ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr); 2863 ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr); 2864 ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 2865 ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr); 2866 for (bc = 0; bc < numBcs; ++bc) { 2867 const PetscInt idx = bcNodes[bc]; 2868 if (idx < n) globalUpdate[idx] = globalRHS[idx]; 2869 } 2870 2871 ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 2872 ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr); 2873 ierr = VecRestoreArray(y, &globalUpdate);CHKERRQ(ierr); 2874 2875 ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr); 2876 ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 2877 PetscFunctionReturn(0); 2878 } 2879 2880 static PetscErrorCode PCReset_PATCH_Linear(PC pc) 2881 { 2882 PC_PATCH *patch = (PC_PATCH *) pc->data; 2883 PetscInt i; 2884 PetscErrorCode ierr; 2885 2886 PetscFunctionBegin; 2887 if (patch->solver) { 2888 for (i = 0; i < patch->npatch; ++i) {ierr = KSPReset((KSP) patch->solver[i]);CHKERRQ(ierr);} 2889 } 2890 PetscFunctionReturn(0); 2891 } 2892 2893 static PetscErrorCode PCReset_PATCH(PC pc) 2894 { 2895 PC_PATCH *patch = (PC_PATCH *) pc->data; 2896 PetscInt i; 2897 PetscErrorCode ierr; 2898 2899 PetscFunctionBegin; 2900 2901 ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr); 2902 ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr); 2903 ierr = PetscSectionDestroy(&patch->pointCounts);CHKERRQ(ierr); 2904 ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr); 2905 ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr); 2906 ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr); 2907 ierr = ISDestroy(&patch->cells);CHKERRQ(ierr); 2908 ierr = ISDestroy(&patch->points);CHKERRQ(ierr); 2909 ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr); 2910 ierr = ISDestroy(&patch->offs);CHKERRQ(ierr); 2911 ierr = PetscSectionDestroy(&patch->patchSection);CHKERRQ(ierr); 2912 ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr); 2913 ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr); 2914 ierr = PetscSectionDestroy(&patch->gtolCountsWithArtificial);CHKERRQ(ierr); 2915 ierr = ISDestroy(&patch->gtolWithArtificial);CHKERRQ(ierr); 2916 ierr = ISDestroy(&patch->dofsWithArtificial);CHKERRQ(ierr); 2917 ierr = ISDestroy(&patch->offsWithArtificial);CHKERRQ(ierr); 2918 ierr = PetscSectionDestroy(&patch->gtolCountsWithAll);CHKERRQ(ierr); 2919 ierr = ISDestroy(&patch->gtolWithAll);CHKERRQ(ierr); 2920 ierr = ISDestroy(&patch->dofsWithAll);CHKERRQ(ierr); 2921 ierr = ISDestroy(&patch->offsWithAll);CHKERRQ(ierr); 2922 ierr = VecDestroy(&patch->cellMats);CHKERRQ(ierr); 2923 ierr = VecDestroy(&patch->intFacetMats);CHKERRQ(ierr); 2924 ierr = ISDestroy(&patch->allCells);CHKERRQ(ierr); 2925 2926 if (patch->dofSection) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);} 2927 ierr = PetscFree(patch->dofSection);CHKERRQ(ierr); 2928 ierr = PetscFree(patch->bs);CHKERRQ(ierr); 2929 ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr); 2930 if (patch->cellNodeMap) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscFree(patch->cellNodeMap[i]);CHKERRQ(ierr);} 2931 ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr); 2932 ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr); 2933 2934 ierr = (*patch->resetsolver)(pc);CHKERRQ(ierr); 2935 2936 if (patch->subspaces_to_exclude) { 2937 PetscHSetIDestroy(&patch->subspaces_to_exclude); 2938 } 2939 2940 ierr = VecDestroy(&patch->localRHS);CHKERRQ(ierr); 2941 ierr = VecDestroy(&patch->localUpdate);CHKERRQ(ierr); 2942 if (patch->patchRHS) { 2943 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHS[i]);CHKERRQ(ierr);} 2944 ierr = PetscFree(patch->patchRHS);CHKERRQ(ierr); 2945 } 2946 if (patch->patchUpdate) { 2947 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchUpdate[i]);CHKERRQ(ierr);} 2948 ierr = PetscFree(patch->patchUpdate);CHKERRQ(ierr); 2949 } 2950 ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr); 2951 if (patch->patch_dof_weights) { 2952 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);} 2953 ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr); 2954 } 2955 if (patch->mat) { 2956 for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);} 2957 ierr = PetscFree(patch->mat);CHKERRQ(ierr); 2958 } 2959 if (patch->matWithArtificial) { 2960 for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->matWithArtificial[i]);CHKERRQ(ierr);} 2961 ierr = PetscFree(patch->matWithArtificial);CHKERRQ(ierr); 2962 } 2963 if (patch->patchRHSWithArtificial) { 2964 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHSWithArtificial[i]);CHKERRQ(ierr);} 2965 ierr = PetscFree(patch->patchRHSWithArtificial);CHKERRQ(ierr); 2966 } 2967 if(patch->dofMappingWithoutToWithArtificial) { 2968 for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->dofMappingWithoutToWithArtificial[i]);CHKERRQ(ierr);} 2969 ierr = PetscFree(patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr); 2970 2971 } 2972 if(patch->dofMappingWithoutToWithAll) { 2973 for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->dofMappingWithoutToWithAll[i]);CHKERRQ(ierr);} 2974 ierr = PetscFree(patch->dofMappingWithoutToWithAll);CHKERRQ(ierr); 2975 2976 } 2977 ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr); 2978 if (patch->userIS) { 2979 for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);} 2980 ierr = PetscFree(patch->userIS);CHKERRQ(ierr); 2981 } 2982 ierr = PetscFree(patch->precomputedTensorLocations);CHKERRQ(ierr); 2983 ierr = PetscFree(patch->precomputedIntFacetTensorLocations);CHKERRQ(ierr); 2984 2985 patch->bs = 0; 2986 patch->cellNodeMap = NULL; 2987 patch->nsubspaces = 0; 2988 ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr); 2989 2990 ierr = PetscViewerDestroy(&patch->viewerSection);CHKERRQ(ierr); 2991 PetscFunctionReturn(0); 2992 } 2993 2994 static PetscErrorCode PCDestroy_PATCH_Linear(PC pc) 2995 { 2996 PC_PATCH *patch = (PC_PATCH *) pc->data; 2997 PetscInt i; 2998 PetscErrorCode ierr; 2999 3000 PetscFunctionBegin; 3001 if (patch->solver) { 3002 for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy((KSP *) &patch->solver[i]);CHKERRQ(ierr);} 3003 ierr = PetscFree(patch->solver);CHKERRQ(ierr); 3004 } 3005 PetscFunctionReturn(0); 3006 } 3007 3008 static PetscErrorCode PCDestroy_PATCH(PC pc) 3009 { 3010 PC_PATCH *patch = (PC_PATCH *) pc->data; 3011 PetscErrorCode ierr; 3012 3013 PetscFunctionBegin; 3014 ierr = PCReset_PATCH(pc);CHKERRQ(ierr); 3015 ierr = (*patch->destroysolver)(pc);CHKERRQ(ierr); 3016 ierr = PetscFree(pc->data);CHKERRQ(ierr); 3017 PetscFunctionReturn(0); 3018 } 3019 3020 static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc) 3021 { 3022 PC_PATCH *patch = (PC_PATCH *) pc->data; 3023 PCPatchConstructType patchConstructionType = PC_PATCH_STAR; 3024 char sub_mat_type[PETSC_MAX_PATH_LEN]; 3025 char option[PETSC_MAX_PATH_LEN]; 3026 const char *prefix; 3027 PetscBool flg, dimflg, codimflg; 3028 MPI_Comm comm; 3029 PetscInt *ifields, nfields, k; 3030 PetscErrorCode ierr; 3031 PCCompositeType loctype = PC_COMPOSITE_ADDITIVE; 3032 3033 PetscFunctionBegin; 3034 ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr); 3035 ierr = PetscObjectGetOptionsPrefix((PetscObject) pc, &prefix);CHKERRQ(ierr); 3036 ierr = PetscOptionsHead(PetscOptionsObject, "Patch solver options");CHKERRQ(ierr); 3037 3038 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_save_operators", patch->classname);CHKERRQ(ierr); 3039 ierr = PetscOptionsBool(option, "Store all patch operators for lifetime of object?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr); 3040 3041 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_precompute_element_tensors", patch->classname); 3042 ierr = PetscOptionsBool(option, "Compute each element tensor only once?", "PCPatchSetPrecomputeElementTensors", patch->precomputeElementTensors, &patch->precomputeElementTensors, &flg);CHKERRQ(ierr); 3043 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_partition_of_unity", patch->classname); 3044 ierr = PetscOptionsBool(option, "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr); 3045 3046 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_local_type", patch->classname);CHKERRQ(ierr); 3047 ierr = PetscOptionsEnum(option,"Type of local solver composition (additive or multiplicative)","PCPatchSetLocalComposition",PCCompositeTypes,(PetscEnum)loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr); 3048 if(flg) { ierr = PCPatchSetLocalComposition(pc, loctype);CHKERRQ(ierr);} 3049 3050 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_dim", patch->classname);CHKERRQ(ierr); 3051 ierr = PetscOptionsInt(option, "What dimension of mesh point to construct patches by? (0 = vertices)", "PCPATCH", patch->dim, &patch->dim, &dimflg);CHKERRQ(ierr); 3052 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_codim", patch->classname);CHKERRQ(ierr); 3053 ierr = PetscOptionsInt(option, "What co-dimension of mesh point to construct patches by? (0 = cells)", "PCPATCH", patch->codim, &patch->codim, &codimflg);CHKERRQ(ierr); 3054 if (dimflg && codimflg) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension"); 3055 3056 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_type", patch->classname);CHKERRQ(ierr); 3057 ierr = PetscOptionsEnum(option, "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr); 3058 if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);} 3059 3060 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_vanka_dim", patch->classname);CHKERRQ(ierr); 3061 ierr = PetscOptionsInt(option, "Topological dimension of entities for Vanka to ignore", "PCPATCH", patch->vankadim, &patch->vankadim, &flg);CHKERRQ(ierr); 3062 3063 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_ignore_dim", patch->classname);CHKERRQ(ierr); 3064 ierr = PetscOptionsInt(option, "Topological dimension of entities for completion to ignore", "PCPATCH", patch->ignoredim, &patch->ignoredim, &flg);CHKERRQ(ierr); 3065 3066 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_pardecomp_overlap", patch->classname);CHKERRQ(ierr); 3067 ierr = PetscOptionsInt(option, "What overlap should we use in construct type pardecomp?", "PCPATCH", patch->pardecomp_overlap, &patch->pardecomp_overlap, &flg);CHKERRQ(ierr); 3068 3069 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_sub_mat_type", patch->classname);CHKERRQ(ierr); 3070 ierr = PetscOptionsFList(option, "Matrix type for patch solves", "PCPatchSetSubMatType", MatList, NULL, sub_mat_type, PETSC_MAX_PATH_LEN, &flg);CHKERRQ(ierr); 3071 if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);} 3072 3073 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_symmetrise_sweep", patch->classname);CHKERRQ(ierr); 3074 ierr = PetscOptionsBool(option, "Go start->end, end->start?", "PCPATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr); 3075 3076 /* If the user has set the number of subspaces, use that for the buffer size, 3077 otherwise use a large number */ 3078 if (patch->nsubspaces <= 0) { 3079 nfields = 128; 3080 } else { 3081 nfields = patch->nsubspaces; 3082 } 3083 ierr = PetscMalloc1(nfields, &ifields);CHKERRQ(ierr); 3084 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_exclude_subspaces", patch->classname);CHKERRQ(ierr); 3085 ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,option,ifields,&nfields,&flg);CHKERRQ(ierr); 3086 if (flg && (patchConstructionType == PC_PATCH_USER)) SETERRQ(comm, PETSC_ERR_ARG_INCOMP, "We cannot support excluding a subspace with user patches because we do not index patches with a mesh point"); 3087 if (flg) { 3088 PetscHSetIClear(patch->subspaces_to_exclude); 3089 for (k = 0; k < nfields; k++) { 3090 PetscHSetIAdd(patch->subspaces_to_exclude, ifields[k]); 3091 } 3092 } 3093 ierr = PetscFree(ifields);CHKERRQ(ierr); 3094 3095 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_patches_view", patch->classname);CHKERRQ(ierr); 3096 ierr = PetscOptionsBool(option, "Print out information during patch construction", "PCPATCH", patch->viewPatches, &patch->viewPatches, &flg);CHKERRQ(ierr); 3097 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_cells_view", patch->classname);CHKERRQ(ierr); 3098 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerCells, &patch->formatCells, &patch->viewCells);CHKERRQ(ierr); 3099 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_interior_facets_view", patch->classname);CHKERRQ(ierr); 3100 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerIntFacets, &patch->formatIntFacets, &patch->viewIntFacets);CHKERRQ(ierr); 3101 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_exterior_facets_view", patch->classname);CHKERRQ(ierr); 3102 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerExtFacets, &patch->formatExtFacets, &patch->viewExtFacets);CHKERRQ(ierr); 3103 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_points_view", patch->classname);CHKERRQ(ierr); 3104 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerPoints, &patch->formatPoints, &patch->viewPoints);CHKERRQ(ierr); 3105 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_section_view", patch->classname);CHKERRQ(ierr); 3106 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerSection, &patch->formatSection, &patch->viewSection);CHKERRQ(ierr); 3107 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_mat_view", patch->classname);CHKERRQ(ierr); 3108 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options, prefix, option, &patch->viewerMatrix, &patch->formatMatrix, &patch->viewMatrix);CHKERRQ(ierr); 3109 ierr = PetscOptionsTail();CHKERRQ(ierr); 3110 patch->optionsSet = PETSC_TRUE; 3111 PetscFunctionReturn(0); 3112 } 3113 3114 static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc) 3115 { 3116 PC_PATCH *patch = (PC_PATCH*) pc->data; 3117 KSPConvergedReason reason; 3118 PetscInt i; 3119 PetscErrorCode ierr; 3120 3121 PetscFunctionBegin; 3122 if (!patch->save_operators) { 3123 /* Can't do this here because the sub KSPs don't have an operator attached yet. */ 3124 PetscFunctionReturn(0); 3125 } 3126 for (i = 0; i < patch->npatch; ++i) { 3127 if (!((KSP) patch->solver[i])->setfromoptionscalled) { 3128 ierr = KSPSetFromOptions((KSP) patch->solver[i]);CHKERRQ(ierr); 3129 } 3130 ierr = KSPSetUp((KSP) patch->solver[i]);CHKERRQ(ierr); 3131 ierr = KSPGetConvergedReason((KSP) patch->solver[i], &reason);CHKERRQ(ierr); 3132 if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 3133 } 3134 PetscFunctionReturn(0); 3135 } 3136 3137 static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer) 3138 { 3139 PC_PATCH *patch = (PC_PATCH *) pc->data; 3140 PetscViewer sviewer; 3141 PetscBool isascii; 3142 PetscMPIInt rank; 3143 PetscErrorCode ierr; 3144 3145 PetscFunctionBegin; 3146 /* TODO Redo tabbing with set tbas in new style */ 3147 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 3148 if (!isascii) PetscFunctionReturn(0); 3149 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr); 3150 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 3151 ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr); 3152 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 3153 ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: multiplicative\n");CHKERRQ(ierr); 3154 } else { 3155 ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr); 3156 } 3157 if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);} 3158 else {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);} 3159 if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);} 3160 else {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);} 3161 if (!patch->precomputeElementTensors) {ierr = PetscViewerASCIIPrintf(viewer, "Not precomputing element tensors (overlapping cells rebuilt in every patch assembly)\n");CHKERRQ(ierr);} 3162 else {ierr = PetscViewerASCIIPrintf(viewer, "Precomputing element tensors (each cell assembled only once)\n");CHKERRQ(ierr);} 3163 if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);} 3164 else {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);} 3165 if (patch->patchconstructop == PCPatchConstruct_Star) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);} 3166 else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);} 3167 else if (patch->patchconstructop == PCPatchConstruct_User) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);} 3168 else {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);} 3169 3170 if (patch->isNonlinear) { 3171 ierr = PetscViewerASCIIPrintf(viewer, "SNES on patches (all same):\n");CHKERRQ(ierr); 3172 } else { 3173 ierr = PetscViewerASCIIPrintf(viewer, "KSP on patches (all same):\n");CHKERRQ(ierr); 3174 } 3175 if (patch->solver) { 3176 ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 3177 if (!rank) { 3178 ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr); 3179 ierr = PetscObjectView(patch->solver[0], sviewer);CHKERRQ(ierr); 3180 ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr); 3181 } 3182 ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 3183 } else { 3184 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 3185 ierr = PetscViewerASCIIPrintf(viewer, "Solver not yet set.\n");CHKERRQ(ierr); 3186 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 3187 } 3188 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 3189 PetscFunctionReturn(0); 3190 } 3191 3192 /*MC 3193 PCPATCH - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping 3194 small block additive preconditioners. Block definition is based on topology from 3195 a DM and equation numbering from a PetscSection. 3196 3197 Options Database Keys: 3198 + -pc_patch_cells_view - Views the process local cell numbers for each patch 3199 . -pc_patch_points_view - Views the process local mesh point numbers for each patch 3200 . -pc_patch_g2l_view - Views the map between global dofs and patch local dofs for each patch 3201 . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary 3202 - -pc_patch_sub_mat_view - Views the matrix associated with each patch 3203 3204 Level: intermediate 3205 3206 .seealso: PCType, PCCreate(), PCSetType() 3207 M*/ 3208 PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc) 3209 { 3210 PC_PATCH *patch; 3211 PetscErrorCode ierr; 3212 3213 PetscFunctionBegin; 3214 ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr); 3215 3216 if (patch->subspaces_to_exclude) { 3217 PetscHSetIDestroy(&patch->subspaces_to_exclude); 3218 } 3219 PetscHSetICreate(&patch->subspaces_to_exclude); 3220 3221 patch->classname = "pc"; 3222 patch->isNonlinear = PETSC_FALSE; 3223 3224 /* Set some defaults */ 3225 patch->combined = PETSC_FALSE; 3226 patch->save_operators = PETSC_TRUE; 3227 patch->local_composition_type = PC_COMPOSITE_ADDITIVE; 3228 patch->precomputeElementTensors = PETSC_FALSE; 3229 patch->partition_of_unity = PETSC_FALSE; 3230 patch->codim = -1; 3231 patch->dim = -1; 3232 patch->vankadim = -1; 3233 patch->ignoredim = -1; 3234 patch->pardecomp_overlap = 0; 3235 patch->patchconstructop = PCPatchConstruct_Star; 3236 patch->symmetrise_sweep = PETSC_FALSE; 3237 patch->npatch = 0; 3238 patch->userIS = NULL; 3239 patch->optionsSet = PETSC_FALSE; 3240 patch->iterationSet = NULL; 3241 patch->user_patches = PETSC_FALSE; 3242 ierr = PetscStrallocpy(MATDENSE, (char **) &patch->sub_mat_type);CHKERRQ(ierr); 3243 patch->viewPatches = PETSC_FALSE; 3244 patch->viewCells = PETSC_FALSE; 3245 patch->viewPoints = PETSC_FALSE; 3246 patch->viewSection = PETSC_FALSE; 3247 patch->viewMatrix = PETSC_FALSE; 3248 patch->setupsolver = PCSetUp_PATCH_Linear; 3249 patch->applysolver = PCApply_PATCH_Linear; 3250 patch->resetsolver = PCReset_PATCH_Linear; 3251 patch->destroysolver = PCDestroy_PATCH_Linear; 3252 patch->updatemultiplicative = PCUpdateMultiplicative_PATCH_Linear; 3253 patch->dofMappingWithoutToWithArtificial = NULL; 3254 patch->dofMappingWithoutToWithAll = NULL; 3255 3256 pc->data = (void *) patch; 3257 pc->ops->apply = PCApply_PATCH; 3258 pc->ops->applytranspose = 0; /* PCApplyTranspose_PATCH; */ 3259 pc->ops->setup = PCSetUp_PATCH; 3260 pc->ops->reset = PCReset_PATCH; 3261 pc->ops->destroy = PCDestroy_PATCH; 3262 pc->ops->setfromoptions = PCSetFromOptions_PATCH; 3263 pc->ops->setuponblocks = PCSetUpOnBlocks_PATCH; 3264 pc->ops->view = PCView_PATCH; 3265 pc->ops->applyrichardson = 0; 3266 3267 PetscFunctionReturn(0); 3268 } 3269