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