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 /* The user's already set the patches in patch->userIS. Build the hash tables */ 80 static PetscErrorCode PCPatchConstruct_User(void *vpatch, DM dm, PetscInt point, PetscHSetI ht) 81 { 82 PC_PATCH *patch = (PC_PATCH *) vpatch; 83 IS patchis = patch->userIS[point]; 84 PetscInt n; 85 const PetscInt *patchdata; 86 PetscInt pStart, pEnd, i; 87 PetscErrorCode ierr; 88 89 PetscFunctionBegin; 90 ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 91 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 92 ierr = ISGetLocalSize(patchis, &n);CHKERRQ(ierr); 93 ierr = ISGetIndices(patchis, &patchdata);CHKERRQ(ierr); 94 for (i = 0; i < n; ++i) { 95 const PetscInt ownedpoint = patchdata[i]; 96 97 if (ownedpoint < pStart || ownedpoint >= pEnd) { 98 SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D was not in [%D, %D)", ownedpoint, pStart, pEnd); 99 } 100 ierr = PetscHSetIAdd(ht, ownedpoint);CHKERRQ(ierr); 101 } 102 ierr = ISRestoreIndices(patchis, &patchdata);CHKERRQ(ierr); 103 PetscFunctionReturn(0); 104 } 105 106 static PetscErrorCode PCPatchCreateDefaultSF_Private(PC pc, PetscInt n, const PetscSF *sf, const PetscInt *bs) 107 { 108 PC_PATCH *patch = (PC_PATCH *) pc->data; 109 PetscInt i; 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 if (n == 1 && bs[0] == 1) { 114 patch->defaultSF = sf[0]; 115 ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr); 116 } else { 117 PetscInt allRoots = 0, allLeaves = 0; 118 PetscInt leafOffset = 0; 119 PetscInt *ilocal = NULL; 120 PetscSFNode *iremote = NULL; 121 PetscInt *remoteOffsets = NULL; 122 PetscInt index = 0; 123 PetscHMapI rankToIndex; 124 PetscInt numRanks = 0; 125 PetscSFNode *remote = NULL; 126 PetscSF rankSF; 127 PetscInt *ranks = NULL; 128 PetscInt *offsets = NULL; 129 MPI_Datatype contig; 130 PetscHSetI ranksUniq; 131 132 /* First figure out how many dofs there are in the concatenated numbering. 133 * allRoots: number of owned global dofs; 134 * allLeaves: number of visible dofs (global + ghosted). 135 */ 136 for (i = 0; i < n; ++i) { 137 PetscInt nroots, nleaves; 138 139 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr); 140 allRoots += nroots * bs[i]; 141 allLeaves += nleaves * bs[i]; 142 } 143 ierr = PetscMalloc1(allLeaves, &ilocal);CHKERRQ(ierr); 144 ierr = PetscMalloc1(allLeaves, &iremote);CHKERRQ(ierr); 145 /* Now build an SF that just contains process connectivity. */ 146 ierr = PetscHSetICreate(&ranksUniq);CHKERRQ(ierr); 147 for (i = 0; i < n; ++i) { 148 const PetscMPIInt *ranks = NULL; 149 PetscInt nranks, j; 150 151 ierr = PetscSFSetUp(sf[i]);CHKERRQ(ierr); 152 ierr = PetscSFGetRanks(sf[i], &nranks, &ranks, NULL, NULL, NULL);CHKERRQ(ierr); 153 /* These are all the ranks who communicate with me. */ 154 for (j = 0; j < nranks; ++j) { 155 ierr = PetscHSetIAdd(ranksUniq, (PetscInt) ranks[j]);CHKERRQ(ierr); 156 } 157 } 158 ierr = PetscHSetIGetSize(ranksUniq, &numRanks);CHKERRQ(ierr); 159 ierr = PetscMalloc1(numRanks, &remote);CHKERRQ(ierr); 160 ierr = PetscMalloc1(numRanks, &ranks);CHKERRQ(ierr); 161 ierr = PetscHSetIGetElems(ranksUniq, &index, ranks);CHKERRQ(ierr); 162 163 ierr = PetscHMapICreate(&rankToIndex);CHKERRQ(ierr); 164 for (i = 0; i < numRanks; ++i) { 165 remote[i].rank = ranks[i]; 166 remote[i].index = 0; 167 ierr = PetscHMapISet(rankToIndex, ranks[i], i);CHKERRQ(ierr); 168 } 169 ierr = PetscFree(ranks);CHKERRQ(ierr); 170 ierr = PetscHSetIDestroy(&ranksUniq);CHKERRQ(ierr); 171 ierr = PetscSFCreate(PetscObjectComm((PetscObject) pc), &rankSF);CHKERRQ(ierr); 172 ierr = PetscSFSetGraph(rankSF, 1, numRanks, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr); 173 ierr = PetscSFSetUp(rankSF);CHKERRQ(ierr); 174 /* OK, use it to communicate the root offset on the remote 175 * processes for each subspace. */ 176 ierr = PetscMalloc1(n, &offsets);CHKERRQ(ierr); 177 ierr = PetscMalloc1(n*numRanks, &remoteOffsets);CHKERRQ(ierr); 178 179 offsets[0] = 0; 180 for (i = 1; i < n; ++i) { 181 PetscInt nroots; 182 183 ierr = PetscSFGetGraph(sf[i-1], &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 184 offsets[i] = offsets[i-1] + nroots*bs[i-1]; 185 } 186 /* Offsets are the offsets on the current process of the 187 * global dof numbering for the subspaces. */ 188 ierr = MPI_Type_contiguous(n, MPIU_INT, &contig);CHKERRQ(ierr); 189 ierr = MPI_Type_commit(&contig);CHKERRQ(ierr); 190 191 ierr = PetscSFBcastBegin(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 192 ierr = PetscSFBcastEnd(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr); 193 ierr = MPI_Type_free(&contig);CHKERRQ(ierr); 194 ierr = PetscFree(offsets);CHKERRQ(ierr); 195 ierr = PetscSFDestroy(&rankSF);CHKERRQ(ierr); 196 /* Now remoteOffsets contains the offsets on the remote 197 * processes who communicate with me. So now we can 198 * concatenate the list of SFs into a single one. */ 199 index = 0; 200 for (i = 0; i < n; ++i) { 201 const PetscSFNode *remote = NULL; 202 const PetscInt *local = NULL; 203 PetscInt nroots, nleaves, j; 204 205 ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, &local, &remote);CHKERRQ(ierr); 206 for (j = 0; j < nleaves; ++j) { 207 PetscInt rank = remote[j].rank; 208 PetscInt idx, rootOffset, k; 209 210 ierr = PetscHMapIGet(rankToIndex, rank, &idx);CHKERRQ(ierr); 211 if (idx == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Didn't find rank, huh?"); 212 /* Offset on given rank for ith subspace */ 213 rootOffset = remoteOffsets[n*idx + i]; 214 for (k = 0; k < bs[i]; ++k) { 215 ilocal[index] = (local ? local[j] : j)*bs[i] + k + leafOffset; 216 iremote[index].rank = remote[j].rank; 217 iremote[index].index = remote[j].index*bs[i] + k + rootOffset; 218 ++index; 219 } 220 } 221 leafOffset += nleaves * bs[i]; 222 } 223 ierr = PetscHMapIDestroy(&rankToIndex);CHKERRQ(ierr); 224 ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 225 ierr = PetscSFCreate(PetscObjectComm((PetscObject)pc), &patch->defaultSF);CHKERRQ(ierr); 226 ierr = PetscSFSetGraph(patch->defaultSF, allRoots, allLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 227 } 228 PetscFunctionReturn(0); 229 } 230 231 /* TODO: Docs */ 232 PetscErrorCode PCPatchSetIgnoreDim(PC pc, PetscInt dim) 233 { 234 PC_PATCH *patch = (PC_PATCH *) pc->data; 235 PetscFunctionBegin; 236 patch->ignoredim = dim; 237 PetscFunctionReturn(0); 238 } 239 240 /* TODO: Docs */ 241 PetscErrorCode PCPatchGetIgnoreDim(PC pc, PetscInt *dim) 242 { 243 PC_PATCH *patch = (PC_PATCH *) pc->data; 244 PetscFunctionBegin; 245 *dim = patch->ignoredim; 246 PetscFunctionReturn(0); 247 } 248 249 /* TODO: Docs */ 250 PetscErrorCode PCPatchSetSaveOperators(PC pc, PetscBool flg) 251 { 252 PC_PATCH *patch = (PC_PATCH *) pc->data; 253 PetscFunctionBegin; 254 patch->save_operators = flg; 255 PetscFunctionReturn(0); 256 } 257 258 /* TODO: Docs */ 259 PetscErrorCode PCPatchGetSaveOperators(PC pc, PetscBool *flg) 260 { 261 PC_PATCH *patch = (PC_PATCH *) pc->data; 262 PetscFunctionBegin; 263 *flg = patch->save_operators; 264 PetscFunctionReturn(0); 265 } 266 267 /* TODO: Docs */ 268 PetscErrorCode PCPatchSetPartitionOfUnity(PC pc, PetscBool flg) 269 { 270 PC_PATCH *patch = (PC_PATCH *) pc->data; 271 PetscFunctionBegin; 272 patch->partition_of_unity = flg; 273 PetscFunctionReturn(0); 274 } 275 276 /* TODO: Docs */ 277 PetscErrorCode PCPatchGetPartitionOfUnity(PC pc, PetscBool *flg) 278 { 279 PC_PATCH *patch = (PC_PATCH *) pc->data; 280 PetscFunctionBegin; 281 *flg = patch->partition_of_unity; 282 PetscFunctionReturn(0); 283 } 284 285 /* TODO: Docs */ 286 PetscErrorCode PCPatchSetLocalComposition(PC pc, PCCompositeType type) 287 { 288 PC_PATCH *patch = (PC_PATCH *) pc->data; 289 PetscFunctionBegin; 290 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"); 291 patch->local_composition_type = type; 292 PetscFunctionReturn(0); 293 } 294 295 /* TODO: Docs */ 296 PetscErrorCode PCPatchGetLocalComposition(PC pc, PCCompositeType *type) 297 { 298 PC_PATCH *patch = (PC_PATCH *) pc->data; 299 PetscFunctionBegin; 300 *type = patch->local_composition_type; 301 PetscFunctionReturn(0); 302 } 303 304 /* TODO: Docs */ 305 PetscErrorCode PCPatchSetSubMatType(PC pc, MatType sub_mat_type) 306 { 307 PC_PATCH *patch = (PC_PATCH *) pc->data; 308 PetscErrorCode ierr; 309 310 PetscFunctionBegin; 311 if (patch->sub_mat_type) {ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);} 312 ierr = PetscStrallocpy(sub_mat_type, (char **) &patch->sub_mat_type);CHKERRQ(ierr); 313 PetscFunctionReturn(0); 314 } 315 316 /* TODO: Docs */ 317 PetscErrorCode PCPatchGetSubMatType(PC pc, MatType *sub_mat_type) 318 { 319 PC_PATCH *patch = (PC_PATCH *) pc->data; 320 PetscFunctionBegin; 321 *sub_mat_type = patch->sub_mat_type; 322 PetscFunctionReturn(0); 323 } 324 325 /* TODO: Docs */ 326 PetscErrorCode PCPatchSetCellNumbering(PC pc, PetscSection cellNumbering) 327 { 328 PC_PATCH *patch = (PC_PATCH *) pc->data; 329 PetscErrorCode ierr; 330 331 PetscFunctionBegin; 332 patch->cellNumbering = cellNumbering; 333 ierr = PetscObjectReference((PetscObject) cellNumbering);CHKERRQ(ierr); 334 PetscFunctionReturn(0); 335 } 336 337 /* TODO: Docs */ 338 PetscErrorCode PCPatchGetCellNumbering(PC pc, PetscSection *cellNumbering) 339 { 340 PC_PATCH *patch = (PC_PATCH *) pc->data; 341 PetscFunctionBegin; 342 *cellNumbering = patch->cellNumbering; 343 PetscFunctionReturn(0); 344 } 345 346 /* TODO: Docs */ 347 PetscErrorCode PCPatchSetConstructType(PC pc, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx) 348 { 349 PC_PATCH *patch = (PC_PATCH *) pc->data; 350 351 PetscFunctionBegin; 352 patch->ctype = ctype; 353 switch (ctype) { 354 case PC_PATCH_STAR: 355 patch->user_patches = PETSC_FALSE; 356 patch->patchconstructop = PCPatchConstruct_Star; 357 break; 358 case PC_PATCH_VANKA: 359 patch->user_patches = PETSC_FALSE; 360 patch->patchconstructop = PCPatchConstruct_Vanka; 361 break; 362 case PC_PATCH_USER: 363 case PC_PATCH_PYTHON: 364 patch->user_patches = PETSC_TRUE; 365 patch->patchconstructop = PCPatchConstruct_User; 366 if (func) { 367 patch->userpatchconstructionop = func; 368 patch->userpatchconstructctx = ctx; 369 } 370 break; 371 default: 372 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 373 } 374 PetscFunctionReturn(0); 375 } 376 377 /* TODO: Docs */ 378 PetscErrorCode PCPatchGetConstructType(PC pc, PCPatchConstructType *ctype, PetscErrorCode (**func)(PC, PetscInt *, IS **, IS *, void *), void **ctx) 379 { 380 PC_PATCH *patch = (PC_PATCH *) pc->data; 381 382 PetscFunctionBegin; 383 *ctype = patch->ctype; 384 switch (patch->ctype) { 385 case PC_PATCH_STAR: 386 case PC_PATCH_VANKA: 387 break; 388 case PC_PATCH_USER: 389 case PC_PATCH_PYTHON: 390 *func = patch->userpatchconstructionop; 391 *ctx = patch->userpatchconstructctx; 392 break; 393 default: 394 SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype); 395 } 396 PetscFunctionReturn(0); 397 } 398 399 /* TODO: Docs */ 400 PetscErrorCode PCPatchSetDiscretisationInfo(PC pc, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, 401 const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes) 402 { 403 PC_PATCH *patch = (PC_PATCH *) pc->data; 404 DM dm; 405 PetscSF *sfs; 406 PetscInt cStart, cEnd, i, j; 407 PetscErrorCode ierr; 408 409 PetscFunctionBegin; 410 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 411 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 412 ierr = PetscMalloc1(nsubspaces, &sfs);CHKERRQ(ierr); 413 ierr = PetscMalloc1(nsubspaces, &patch->dofSection);CHKERRQ(ierr); 414 ierr = PetscMalloc1(nsubspaces, &patch->bs);CHKERRQ(ierr); 415 ierr = PetscMalloc1(nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr); 416 ierr = PetscMalloc1(nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr); 417 ierr = PetscMalloc1(nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr); 418 419 patch->nsubspaces = nsubspaces; 420 patch->totalDofsPerCell = 0; 421 for (i = 0; i < nsubspaces; ++i) { 422 ierr = DMGetDefaultSection(dms[i], &patch->dofSection[i]);CHKERRQ(ierr); 423 ierr = PetscObjectReference((PetscObject) patch->dofSection[i]);CHKERRQ(ierr); 424 ierr = DMGetDefaultSF(dms[i], &sfs[i]);CHKERRQ(ierr); 425 patch->bs[i] = bs[i]; 426 patch->nodesPerCell[i] = nodesPerCell[i]; 427 patch->totalDofsPerCell += nodesPerCell[i]*bs[i]; 428 ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr); 429 for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j]; 430 patch->subspaceOffsets[i] = subspaceOffsets[i]; 431 } 432 ierr = PCPatchCreateDefaultSF_Private(pc, nsubspaces, sfs, patch->bs);CHKERRQ(ierr); 433 ierr = PetscFree(sfs);CHKERRQ(ierr); 434 435 patch->subspaceOffsets[nsubspaces] = subspaceOffsets[nsubspaces]; 436 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr); 437 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr); 438 PetscFunctionReturn(0); 439 } 440 441 /* TODO: Docs */ 442 PetscErrorCode PCPatchSetDiscretisationInfoCombined(PC pc, DM dm, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes) 443 { 444 PC_PATCH *patch = (PC_PATCH *) pc->data; 445 PetscInt cStart, cEnd, i, j; 446 PetscErrorCode ierr; 447 448 PetscFunctionBegin; 449 patch->combined = PETSC_TRUE; 450 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 451 ierr = DMGetNumFields(dm, &patch->nsubspaces);CHKERRQ(ierr); 452 ierr = PetscCalloc1(patch->nsubspaces, &patch->dofSection);CHKERRQ(ierr); 453 ierr = PetscMalloc1(patch->nsubspaces, &patch->bs);CHKERRQ(ierr); 454 ierr = PetscMalloc1(patch->nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr); 455 ierr = PetscMalloc1(patch->nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr); 456 ierr = PetscCalloc1(patch->nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr); 457 ierr = DMGetDefaultSection(dm, &patch->dofSection[0]);CHKERRQ(ierr); 458 ierr = PetscObjectReference((PetscObject) patch->dofSection[0]);CHKERRQ(ierr); 459 ierr = PetscSectionGetStorageSize(patch->dofSection[0], &patch->subspaceOffsets[patch->nsubspaces]);CHKERRQ(ierr); 460 patch->totalDofsPerCell = 0; 461 for (i = 0; i < patch->nsubspaces; ++i) { 462 patch->bs[i] = 1; 463 patch->nodesPerCell[i] = nodesPerCell[i]; 464 patch->totalDofsPerCell += nodesPerCell[i]; 465 ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr); 466 for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j]; 467 } 468 ierr = DMGetDefaultSF(dm, &patch->defaultSF);CHKERRQ(ierr); 469 ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr); 470 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr); 471 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr); 472 PetscFunctionReturn(0); 473 } 474 475 /*@C 476 477 PCPatchSetComputeFunction - Set the callback used to compute patch residuals 478 479 Input Parameters: 480 + pc - The PC 481 . func - The callback 482 - ctx - The user context 483 484 Level: advanced 485 486 Note: 487 The callback has signature: 488 + usercomputef(pc, point, x, f, cellIS, n, u, ctx) 489 + pc - The PC 490 + point - The point 491 + x - The input solution (not used in linear problems) 492 + f - The patch residual vector 493 + cellIS - An array of the cell numbers 494 + n - The size of g2l 495 + g2l - The global to local dof translation table 496 + ctx - The user context 497 and can assume that the matrix entries have been set to zero before the call. 498 499 .seealso: PCPatchSetComputeOperator(), PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo() 500 @*/ 501 PetscErrorCode PCPatchSetComputeFunction(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 502 { 503 PC_PATCH *patch = (PC_PATCH *) pc->data; 504 505 PetscFunctionBegin; 506 patch->usercomputef = func; 507 patch->usercomputefctx = ctx; 508 PetscFunctionReturn(0); 509 } 510 511 /*@C 512 513 PCPatchSetComputeOperator - Set the callback used to compute patch matrices 514 515 Input Parameters: 516 + pc - The PC 517 . func - The callback 518 - ctx - The user context 519 520 Level: advanced 521 522 Note: 523 The callback has signature: 524 + usercomputeop(pc, point, x, mat, cellIS, n, u, ctx) 525 + pc - The PC 526 + point - The point 527 + x - The input solution (not used in linear problems) 528 + mat - The patch matrix 529 + cellIS - An array of the cell numbers 530 + n - The size of g2l 531 + g2l - The global to local dof translation table 532 + ctx - The user context 533 and can assume that the matrix entries have been set to zero before the call. 534 535 .seealso: PCPatchGetComputeOperator(), PCPatchSetComputeFunction(), PCPatchSetDiscretisationInfo() 536 @*/ 537 PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *ctx) 538 { 539 PC_PATCH *patch = (PC_PATCH *) pc->data; 540 541 PetscFunctionBegin; 542 patch->usercomputeop = func; 543 patch->usercomputeopctx = ctx; 544 PetscFunctionReturn(0); 545 } 546 547 /* On entry, ht contains the topological entities whose dofs we are responsible for solving for; 548 on exit, cht contains all the topological entities we need to compute their residuals. 549 In full generality this should incorporate knowledge of the sparsity pattern of the matrix; 550 here we assume a standard FE sparsity pattern.*/ 551 /* TODO: Use DMPlexGetAdjacency() */ 552 static PetscErrorCode PCPatchCompleteCellPatch(PC pc, PetscHSetI ht, PetscHSetI cht) 553 { 554 DM dm; 555 PetscHashIter hi; 556 PetscInt point; 557 PetscInt *star = NULL, *closure = NULL; 558 PetscInt ignoredim, iStart = 0, iEnd = -1, starSize, closureSize, si, ci; 559 PetscErrorCode ierr; 560 561 PetscFunctionBegin; 562 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 563 ierr = PCPatchGetIgnoreDim(pc, &ignoredim);CHKERRQ(ierr); 564 if (ignoredim >= 0) {ierr = DMPlexGetDepthStratum(dm, ignoredim, &iStart, &iEnd);CHKERRQ(ierr);} 565 ierr = PetscHSetIClear(cht);CHKERRQ(ierr); 566 PetscHashIterBegin(ht, hi); 567 while (!PetscHashIterAtEnd(ht, hi)) { 568 569 PetscHashIterGetKey(ht, hi, point); 570 PetscHashIterNext(ht, hi); 571 572 /* Loop over all the cells that this point connects to */ 573 ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr); 574 for (si = 0; si < starSize*2; si += 2) { 575 const PetscInt ownedpoint = star[si]; 576 /* TODO Check for point in cht before running through closure again */ 577 /* now loop over all entities in the closure of that cell */ 578 ierr = DMPlexGetTransitiveClosure(dm, ownedpoint, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 579 for (ci = 0; ci < closureSize*2; ci += 2) { 580 const PetscInt seenpoint = closure[ci]; 581 if (ignoredim >= 0 && seenpoint >= iStart && seenpoint < iEnd) continue; 582 ierr = PetscHSetIAdd(cht, seenpoint);CHKERRQ(ierr); 583 } 584 } 585 } 586 ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr); 587 ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_FALSE, NULL, &star);CHKERRQ(ierr); 588 PetscFunctionReturn(0); 589 } 590 591 static PetscErrorCode PCPatchGetGlobalDofs(PC pc, PetscSection dofSection[], PetscInt f, PetscBool combined, PetscInt p, PetscInt *dof, PetscInt *off) 592 { 593 PetscErrorCode ierr; 594 595 PetscFunctionBegin; 596 if (combined) { 597 if (f < 0) { 598 if (dof) {ierr = PetscSectionGetDof(dofSection[0], p, dof);CHKERRQ(ierr);} 599 if (off) {ierr = PetscSectionGetOffset(dofSection[0], p, off);CHKERRQ(ierr);} 600 } else { 601 if (dof) {ierr = PetscSectionGetFieldDof(dofSection[0], p, f, dof);CHKERRQ(ierr);} 602 if (off) {ierr = PetscSectionGetFieldOffset(dofSection[0], p, f, off);CHKERRQ(ierr);} 603 } 604 } else { 605 if (f < 0) { 606 PC_PATCH *patch = (PC_PATCH *) pc->data; 607 PetscInt fdof, g; 608 609 if (dof) { 610 *dof = 0; 611 for (g = 0; g < patch->nsubspaces; ++g) { 612 ierr = PetscSectionGetDof(dofSection[g], p, &fdof);CHKERRQ(ierr); 613 *dof += fdof; 614 } 615 } 616 if (off) { 617 *off = 0; 618 for (g = 0; g < patch->nsubspaces; ++g) { 619 ierr = PetscSectionGetOffset(dofSection[g], p, &fdof);CHKERRQ(ierr); 620 *off += fdof; 621 } 622 } 623 } else { 624 if (dof) {ierr = PetscSectionGetDof(dofSection[f], p, dof);CHKERRQ(ierr);} 625 if (off) {ierr = PetscSectionGetOffset(dofSection[f], p, off);CHKERRQ(ierr);} 626 } 627 } 628 PetscFunctionReturn(0); 629 } 630 631 /* Given a hash table with a set of topological entities (pts), compute the degrees of 632 freedom in global concatenated numbering on those entities. 633 For Vanka smoothing, this needs to do something special: ignore dofs of the 634 constraint subspace on entities that aren't the base entity we're building the patch 635 around. */ 636 static PetscErrorCode PCPatchGetPointDofs(PC pc, PetscHSetI pts, PetscHSetI dofs, PetscInt base, PetscHSetI* subspaces_to_exclude) 637 { 638 PC_PATCH *patch = (PC_PATCH *) pc->data; 639 PetscHashIter hi; 640 PetscInt ldof, loff; 641 PetscInt k, p; 642 PetscErrorCode ierr; 643 644 PetscFunctionBegin; 645 ierr = PetscHSetIClear(dofs);CHKERRQ(ierr); 646 for (k = 0; k < patch->nsubspaces; ++k) { 647 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 648 PetscInt bs = patch->bs[k]; 649 PetscInt j, l; 650 651 if (subspaces_to_exclude != NULL) { 652 PetscBool should_exclude_k = PETSC_FALSE; 653 PetscHSetIHas(*subspaces_to_exclude, k, &should_exclude_k); 654 if (should_exclude_k) { 655 /* only get this subspace dofs at the base entity, not any others */ 656 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, base, &ldof, &loff);CHKERRQ(ierr); 657 if (0 == ldof) continue; 658 for (j = loff; j < ldof + loff; ++j) { 659 for (l = 0; l < bs; ++l) { 660 PetscInt dof = bs*j + l + subspaceOffset; 661 ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr); 662 } 663 } 664 continue; /* skip the other dofs of this subspace */ 665 } 666 } 667 668 PetscHashIterBegin(pts, hi); 669 while (!PetscHashIterAtEnd(pts, hi)) { 670 PetscHashIterGetKey(pts, hi, p); 671 PetscHashIterNext(pts, hi); 672 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, p, &ldof, &loff);CHKERRQ(ierr); 673 if (0 == ldof) continue; 674 for (j = loff; j < ldof + loff; ++j) { 675 for (l = 0; l < bs; ++l) { 676 PetscInt dof = bs*j + l + subspaceOffset; 677 ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr); 678 } 679 } 680 } 681 } 682 PetscFunctionReturn(0); 683 } 684 685 /* Given two hash tables A and B, compute the keys in B that are not in A, and put them in C */ 686 static PetscErrorCode PCPatchComputeSetDifference_Private(PetscHSetI A, PetscHSetI B, PetscHSetI C) 687 { 688 PetscHashIter hi; 689 PetscInt key; 690 PetscBool flg; 691 PetscErrorCode ierr; 692 693 PetscFunctionBegin; 694 ierr = PetscHSetIClear(C);CHKERRQ(ierr); 695 PetscHashIterBegin(B, hi); 696 while (!PetscHashIterAtEnd(B, hi)) { 697 PetscHashIterGetKey(B, hi, key); 698 PetscHashIterNext(B, hi); 699 ierr = PetscHSetIHas(A, key, &flg);CHKERRQ(ierr); 700 if (!flg) {ierr = PetscHSetIAdd(C, key);CHKERRQ(ierr);} 701 } 702 PetscFunctionReturn(0); 703 } 704 705 /* 706 * PCPatchCreateCellPatches - create patches. 707 * 708 * Input Parameters: 709 * + dm - The DMPlex object defining the mesh 710 * 711 * Output Parameters: 712 * + cellCounts - Section with counts of cells around each vertex 713 * . cells - IS of the cell point indices of cells in each patch 714 * . pointCounts - Section with counts of cells around each vertex 715 * - point - IS of the cell point indices of cells in each patch 716 */ 717 static PetscErrorCode PCPatchCreateCellPatches(PC pc) 718 { 719 PC_PATCH *patch = (PC_PATCH *) pc->data; 720 DMLabel ghost = NULL; 721 DM dm, plex; 722 PetscHSetI ht, cht; 723 PetscSection cellCounts, pointCounts; 724 PetscInt *cellsArray, *pointsArray; 725 PetscInt numCells, numPoints; 726 const PetscInt *leaves; 727 PetscInt nleaves, pStart, pEnd, cStart, cEnd, vStart, vEnd, v; 728 PetscBool isFiredrake; 729 PetscErrorCode ierr; 730 731 PetscFunctionBegin; 732 /* Used to keep track of the cells in the patch. */ 733 ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 734 ierr = PetscHSetICreate(&cht);CHKERRQ(ierr); 735 736 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 737 if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch PC\n"); 738 ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 739 ierr = DMPlexGetChart(plex, &pStart, &pEnd);CHKERRQ(ierr); 740 ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr); 741 742 if (patch->user_patches) { 743 ierr = patch->userpatchconstructionop(pc, &patch->npatch, &patch->userIS, &patch->iterationSet, patch->userpatchconstructctx);CHKERRQ(ierr); 744 vStart = 0; vEnd = patch->npatch; 745 } else if (patch->codim < 0) { 746 if (patch->dim < 0) {ierr = DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd);CHKERRQ(ierr);} 747 else {ierr = DMPlexGetDepthStratum(plex, patch->dim, &vStart, &vEnd);CHKERRQ(ierr);} 748 } else {ierr = DMPlexGetHeightStratum(plex, patch->codim, &vStart, &vEnd);CHKERRQ(ierr);} 749 patch->npatch = vEnd - vStart; 750 751 /* These labels mark the owned points. We only create patches around points that this process owns. */ 752 ierr = DMHasLabel(dm, "pyop2_ghost", &isFiredrake);CHKERRQ(ierr); 753 if (isFiredrake) { 754 ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr); 755 ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr); 756 } else { 757 PetscSF sf; 758 759 ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 760 ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 761 nleaves = PetscMax(nleaves, 0); 762 } 763 764 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->cellCounts);CHKERRQ(ierr); 765 ierr = PetscObjectSetName((PetscObject) patch->cellCounts, "Patch Cell Layout");CHKERRQ(ierr); 766 cellCounts = patch->cellCounts; 767 ierr = PetscSectionSetChart(cellCounts, vStart, vEnd);CHKERRQ(ierr); 768 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->pointCounts);CHKERRQ(ierr); 769 ierr = PetscObjectSetName((PetscObject) patch->pointCounts, "Patch Point Layout");CHKERRQ(ierr); 770 pointCounts = patch->pointCounts; 771 ierr = PetscSectionSetChart(pointCounts, vStart, vEnd);CHKERRQ(ierr); 772 /* Count cells and points in the patch surrounding each entity */ 773 for (v = vStart; v < vEnd; ++v) { 774 PetscHashIter hi; 775 PetscInt chtSize, loc = -1; 776 PetscBool flg; 777 778 if (!patch->user_patches) { 779 if (ghost) {ierr = DMLabelHasPoint(ghost, v, &flg);CHKERRQ(ierr);} 780 else {ierr = PetscFindInt(v, nleaves, leaves, &loc);CHKERRQ(ierr); flg = loc >=0 ? PETSC_TRUE : PETSC_FALSE;} 781 /* Not an owned entity, don't make a cell patch. */ 782 if (flg) continue; 783 } 784 785 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 786 ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr); 787 ierr = PetscHSetIGetSize(cht, &chtSize);CHKERRQ(ierr); 788 /* empty patch, continue */ 789 if (chtSize == 0) continue; 790 791 /* safe because size(cht) > 0 from above */ 792 PetscHashIterBegin(cht, hi); 793 while (!PetscHashIterAtEnd(cht, hi)) { 794 PetscInt point, pdof; 795 796 PetscHashIterGetKey(cht, hi, point); 797 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr); 798 if (pdof) {ierr = PetscSectionAddDof(pointCounts, v, 1);CHKERRQ(ierr);} 799 if (point >= cStart && point < cEnd) {ierr = PetscSectionAddDof(cellCounts, v, 1);CHKERRQ(ierr);} 800 PetscHashIterNext(cht, hi); 801 } 802 } 803 if (isFiredrake) {ierr = DMLabelDestroyIndex(ghost);CHKERRQ(ierr);} 804 805 ierr = PetscSectionSetUp(cellCounts);CHKERRQ(ierr); 806 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 807 ierr = PetscMalloc1(numCells, &cellsArray);CHKERRQ(ierr); 808 ierr = PetscSectionSetUp(pointCounts);CHKERRQ(ierr); 809 ierr = PetscSectionGetStorageSize(pointCounts, &numPoints);CHKERRQ(ierr); 810 ierr = PetscMalloc1(numPoints, &pointsArray);CHKERRQ(ierr); 811 812 /* Now that we know how much space we need, run through again and actually remember the cells. */ 813 for (v = vStart; v < vEnd; v++ ) { 814 PetscHashIter hi; 815 PetscInt dof, off, cdof, coff, pdof, n = 0, cn = 0; 816 817 ierr = PetscSectionGetDof(pointCounts, v, &dof);CHKERRQ(ierr); 818 ierr = PetscSectionGetOffset(pointCounts, v, &off);CHKERRQ(ierr); 819 ierr = PetscSectionGetDof(cellCounts, v, &cdof);CHKERRQ(ierr); 820 ierr = PetscSectionGetOffset(cellCounts, v, &coff);CHKERRQ(ierr); 821 if (dof <= 0) continue; 822 ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr); 823 ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr); 824 PetscHashIterBegin(cht, hi); 825 while (!PetscHashIterAtEnd(cht, hi)) { 826 PetscInt point; 827 828 PetscHashIterGetKey(cht, hi, point); 829 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr); 830 if (pdof) {pointsArray[off + n++] = point;} 831 if (point >= cStart && point < cEnd) {cellsArray[coff + cn++] = point;} 832 PetscHashIterNext(cht, hi); 833 } 834 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); 835 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); 836 } 837 ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 838 ierr = PetscHSetIDestroy(&cht);CHKERRQ(ierr); 839 ierr = DMDestroy(&plex);CHKERRQ(ierr); 840 841 ierr = ISCreateGeneral(PETSC_COMM_SELF, numCells, cellsArray, PETSC_OWN_POINTER, &patch->cells);CHKERRQ(ierr); 842 ierr = PetscObjectSetName((PetscObject) patch->cells, "Patch Cells");CHKERRQ(ierr); 843 if (patch->viewCells) { 844 ierr = ObjectView((PetscObject) patch->cellCounts, patch->viewerCells, patch->formatCells);CHKERRQ(ierr); 845 ierr = ObjectView((PetscObject) patch->cells, patch->viewerCells, patch->formatCells);CHKERRQ(ierr); 846 } 847 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, pointsArray, PETSC_OWN_POINTER, &patch->points);CHKERRQ(ierr); 848 ierr = PetscObjectSetName((PetscObject) patch->points, "Patch Points");CHKERRQ(ierr); 849 if (patch->viewPoints) { 850 ierr = ObjectView((PetscObject) patch->pointCounts, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr); 851 ierr = ObjectView((PetscObject) patch->points, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr); 852 } 853 PetscFunctionReturn(0); 854 } 855 856 /* 857 * PCPatchCreateCellPatchDiscretisationInfo - Build the dof maps for cell patches 858 * 859 * Input Parameters: 860 * + dm - The DMPlex object defining the mesh 861 * . cellCounts - Section with counts of cells around each vertex 862 * . cells - IS of the cell point indices of cells in each patch 863 * . cellNumbering - Section mapping plex cell points to Firedrake cell indices. 864 * . nodesPerCell - number of nodes per cell. 865 * - cellNodeMap - map from cells to node indices (nodesPerCell * numCells) 866 * 867 * Output Parameters: 868 * + dofs - IS of local dof numbers of each cell in the patch, where local is a patch local numbering 869 * . gtolCounts - Section with counts of dofs per cell patch 870 * - gtol - IS mapping from global dofs to local dofs for each patch. 871 */ 872 static PetscErrorCode PCPatchCreateCellPatchDiscretisationInfo(PC pc) 873 { 874 PC_PATCH *patch = (PC_PATCH *) pc->data; 875 PetscSection cellCounts = patch->cellCounts; 876 PetscSection pointCounts = patch->pointCounts; 877 PetscSection gtolCounts, gtolCountsWithArtificial = NULL; 878 IS cells = patch->cells; 879 IS points = patch->points; 880 PetscSection cellNumbering = patch->cellNumbering; 881 PetscInt Nf = patch->nsubspaces; 882 PetscInt numCells, numPoints; 883 PetscInt numDofs; 884 PetscInt numGlobalDofs, numGlobalDofsWithArtificial; 885 PetscInt totalDofsPerCell = patch->totalDofsPerCell; 886 PetscInt vStart, vEnd, v; 887 const PetscInt *cellsArray, *pointsArray; 888 PetscInt *newCellsArray = NULL; 889 PetscInt *dofsArray = NULL; 890 PetscInt *dofsArrayWithArtificial = NULL; 891 PetscInt *offsArray = NULL; 892 PetscInt *offsArrayWithArtificial = NULL; 893 PetscInt *asmArray = NULL; 894 PetscInt *asmArrayWithArtificial = NULL; 895 PetscInt *globalDofsArray = NULL; 896 PetscInt *globalDofsArrayWithArtificial = NULL; 897 PetscInt globalIndex = 0; 898 PetscInt key = 0; 899 PetscInt asmKey = 0; 900 DM dm = NULL; 901 const PetscInt *bcNodes = NULL; 902 PetscHMapI ht; 903 PetscHMapI htWithArtificial; 904 PetscHSetI globalBcs; 905 PetscInt numBcs; 906 PetscHSetI ownedpts, seenpts, owneddofs, seendofs, artificialbcs; 907 PetscInt pStart, pEnd, p, i; 908 char option[PETSC_MAX_PATH_LEN]; 909 PetscBool isNonlinear; 910 PetscErrorCode ierr; 911 912 PetscFunctionBegin; 913 914 ierr = PCGetDM(pc, &dm); CHKERRQ(ierr); 915 /* dofcounts section is cellcounts section * dofPerCell */ 916 ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr); 917 ierr = PetscSectionGetStorageSize(patch->pointCounts, &numPoints);CHKERRQ(ierr); 918 numDofs = numCells * totalDofsPerCell; 919 ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr); 920 ierr = PetscMalloc1(numPoints*Nf, &offsArray);CHKERRQ(ierr); 921 ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr); 922 ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr); 923 ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr); 924 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr); 925 gtolCounts = patch->gtolCounts; 926 ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr); 927 ierr = PetscObjectSetName((PetscObject) patch->gtolCounts, "Patch Global Index Section");CHKERRQ(ierr); 928 929 ierr = PetscStrncmp(patch->classname, "snes", sizeof("snes"), &isNonlinear); 930 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) 931 { 932 ierr = PetscMalloc1(numPoints*Nf, &offsArrayWithArtificial);CHKERRQ(ierr); 933 ierr = PetscMalloc1(numDofs, &asmArrayWithArtificial);CHKERRQ(ierr); 934 ierr = PetscMalloc1(numDofs, &dofsArrayWithArtificial);CHKERRQ(ierr); 935 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCountsWithArtificial);CHKERRQ(ierr); 936 gtolCountsWithArtificial = patch->gtolCountsWithArtificial; 937 ierr = PetscSectionSetChart(gtolCountsWithArtificial, vStart, vEnd);CHKERRQ(ierr); 938 ierr = PetscObjectSetName((PetscObject) patch->gtolCountsWithArtificial, "Patch Global Index Section Including Artificial BCs");CHKERRQ(ierr); 939 } 940 941 /* Outside the patch loop, get the dofs that are globally-enforced Dirichlet 942 conditions */ 943 ierr = PetscHSetICreate(&globalBcs);CHKERRQ(ierr); 944 ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr); 945 ierr = ISGetSize(patch->ghostBcNodes, &numBcs); CHKERRQ(ierr); 946 for (i = 0; i < numBcs; ++i) { 947 ierr = PetscHSetIAdd(globalBcs, bcNodes[i]);CHKERRQ(ierr); /* these are already in concatenated numbering */ 948 } 949 ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr); 950 ierr = ISDestroy(&patch->ghostBcNodes); CHKERRQ(ierr); /* memory optimisation */ 951 952 /* Hash tables for artificial BC construction */ 953 ierr = PetscHSetICreate(&ownedpts);CHKERRQ(ierr); 954 ierr = PetscHSetICreate(&seenpts);CHKERRQ(ierr); 955 ierr = PetscHSetICreate(&owneddofs);CHKERRQ(ierr); 956 ierr = PetscHSetICreate(&seendofs);CHKERRQ(ierr); 957 ierr = PetscHSetICreate(&artificialbcs);CHKERRQ(ierr); 958 959 ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr); 960 ierr = ISGetIndices(points, &pointsArray);CHKERRQ(ierr); 961 ierr = PetscHMapICreate(&ht);CHKERRQ(ierr); 962 ierr = PetscHMapICreate(&htWithArtificial);CHKERRQ(ierr); 963 for (v = vStart; v < vEnd; ++v) { 964 PetscInt localIndex = 0; 965 PetscInt localIndexWithArtificial = 0; 966 PetscInt dof, off, i, j, k, l; 967 968 ierr = PetscHMapIClear(ht);CHKERRQ(ierr); 969 ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr); 970 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 971 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 972 if (dof <= 0) continue; 973 974 /* Calculate the global numbers of the artificial BC dofs here first */ 975 ierr = patch->patchconstructop((void*)patch, dm, v, ownedpts); CHKERRQ(ierr); 976 ierr = PCPatchCompleteCellPatch(pc, ownedpts, seenpts); CHKERRQ(ierr); 977 ierr = PCPatchGetPointDofs(pc, ownedpts, owneddofs, v, &patch->subspaces_to_exclude); CHKERRQ(ierr); 978 ierr = PCPatchGetPointDofs(pc, seenpts, seendofs, v, NULL); CHKERRQ(ierr); 979 ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs); CHKERRQ(ierr); 980 if (patch->viewPatches) { 981 PetscHSetI globalbcdofs; 982 PetscHashIter hi; 983 MPI_Comm comm = PetscObjectComm((PetscObject)pc); 984 985 ierr = PetscHSetICreate(&globalbcdofs);CHKERRQ(ierr); 986 ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v); CHKERRQ(ierr); 987 PetscHashIterBegin(owneddofs, hi); 988 while (!PetscHashIterAtEnd(owneddofs, hi)) { 989 PetscInt globalDof; 990 991 PetscHashIterGetKey(owneddofs, hi, globalDof); 992 PetscHashIterNext(owneddofs, hi); 993 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 994 } 995 ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr); 996 ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v); CHKERRQ(ierr); 997 PetscHashIterBegin(seendofs, hi); 998 while (!PetscHashIterAtEnd(seendofs, hi)) { 999 PetscInt globalDof; 1000 PetscBool flg; 1001 1002 PetscHashIterGetKey(seendofs, hi, globalDof); 1003 PetscHashIterNext(seendofs, hi); 1004 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 1005 1006 ierr = PetscHSetIHas(globalBcs, globalDof, &flg);CHKERRQ(ierr); 1007 if (flg) {ierr = PetscHSetIAdd(globalbcdofs, globalDof);CHKERRQ(ierr);} 1008 } 1009 ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr); 1010 ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v); CHKERRQ(ierr); 1011 ierr = PetscHSetIGetSize(globalbcdofs, &numBcs);CHKERRQ(ierr); 1012 if (numBcs > 0) { 1013 PetscHashIterBegin(globalbcdofs, hi); 1014 while (!PetscHashIterAtEnd(globalbcdofs, hi)) { 1015 PetscInt globalDof; 1016 PetscHashIterGetKey(globalbcdofs, hi, globalDof); 1017 PetscHashIterNext(globalbcdofs, hi); 1018 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr); 1019 } 1020 } 1021 ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr); 1022 ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr); 1023 ierr = PetscHSetIGetSize(artificialbcs, &numBcs);CHKERRQ(ierr); 1024 if (numBcs > 0) { 1025 PetscHashIterBegin(artificialbcs, hi); 1026 while (!PetscHashIterAtEnd(artificialbcs, hi)) { 1027 PetscInt globalDof; 1028 PetscHashIterGetKey(artificialbcs, hi, globalDof); 1029 PetscHashIterNext(artificialbcs, hi); 1030 ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr); 1031 } 1032 } 1033 ierr = PetscSynchronizedPrintf(comm, "\n\n"); CHKERRQ(ierr); 1034 ierr = PetscHSetIDestroy(&globalbcdofs);CHKERRQ(ierr); 1035 } 1036 for (k = 0; k < patch->nsubspaces; ++k) { 1037 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1038 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1039 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1040 PetscInt bs = patch->bs[k]; 1041 1042 for (i = off; i < off + dof; ++i) { 1043 /* Walk over the cells in this patch. */ 1044 const PetscInt c = cellsArray[i]; 1045 PetscInt cell = c; 1046 1047 /* TODO Change this to an IS */ 1048 if (cellNumbering) { 1049 ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr); 1050 if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c); 1051 ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr); 1052 } 1053 newCellsArray[i] = cell; 1054 for (j = 0; j < nodesPerCell; ++j) { 1055 /* For each global dof, map it into contiguous local storage. */ 1056 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset; 1057 /* finally, loop over block size */ 1058 for (l = 0; l < bs; ++l) { 1059 PetscInt localDof; 1060 PetscBool isGlobalBcDof, isArtificialBcDof; 1061 1062 /* first, check if this is either a globally enforced or locally enforced BC dof */ 1063 ierr = PetscHSetIHas(globalBcs, globalDof + l, &isGlobalBcDof);CHKERRQ(ierr); 1064 ierr = PetscHSetIHas(artificialbcs, globalDof + l, &isArtificialBcDof);CHKERRQ(ierr); 1065 1066 /* if it's either, don't ever give it a local dof number */ 1067 if (isGlobalBcDof || isArtificialBcDof) { 1068 dofsArray[globalIndex] = -1; /* don't use this in assembly in this patch */ 1069 } else { 1070 ierr = PetscHMapIGet(ht, globalDof + l, &localDof);CHKERRQ(ierr); 1071 if (localDof == -1) { 1072 localDof = localIndex++; 1073 ierr = PetscHMapISet(ht, globalDof + l, localDof);CHKERRQ(ierr); 1074 } 1075 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 1076 /* And store. */ 1077 dofsArray[globalIndex] = localDof; 1078 } 1079 1080 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1081 if (isGlobalBcDof) { 1082 dofsArrayWithArtificial[globalIndex] = -1; /* don't use this in assembly in this patch */ 1083 } else { 1084 ierr = PetscHMapIGet(htWithArtificial, globalDof + l, &localDof);CHKERRQ(ierr); 1085 if (localDof == -1) { 1086 localDof = localIndexWithArtificial++; 1087 ierr = PetscHMapISet(htWithArtificial, globalDof + l, localDof);CHKERRQ(ierr); 1088 } 1089 if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs); 1090 /* And store.*/ 1091 dofsArrayWithArtificial[globalIndex] = localDof; 1092 } 1093 } 1094 globalIndex++; 1095 } 1096 } 1097 } 1098 } 1099 /*How many local dofs in this patch? */ 1100 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1101 ierr = PetscHMapIGetSize(htWithArtificial, &dof);CHKERRQ(ierr); 1102 ierr = PetscSectionSetDof(gtolCountsWithArtificial, v, dof);CHKERRQ(ierr); 1103 } 1104 ierr = PetscHMapIGetSize(ht, &dof);CHKERRQ(ierr); 1105 ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr); 1106 } 1107 if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex); 1108 ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr); 1109 ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr); 1110 ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr); 1111 1112 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1113 ierr = PetscSectionSetUp(gtolCountsWithArtificial);CHKERRQ(ierr); 1114 ierr = PetscSectionGetStorageSize(gtolCountsWithArtificial, &numGlobalDofsWithArtificial);CHKERRQ(ierr); 1115 ierr = PetscMalloc1(numGlobalDofsWithArtificial, &globalDofsArrayWithArtificial);CHKERRQ(ierr); 1116 } 1117 /* Now populate the global to local map. This could be merged into the above loop if we were willing to deal with reallocs. */ 1118 for (v = vStart; v < vEnd; ++v) { 1119 PetscHashIter hi; 1120 PetscInt dof, off, Np, ooff, i, j, k, l; 1121 1122 ierr = PetscHMapIClear(ht);CHKERRQ(ierr); 1123 ierr = PetscHMapIClear(htWithArtificial);CHKERRQ(ierr); 1124 ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr); 1125 ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr); 1126 ierr = PetscSectionGetDof(pointCounts, v, &Np);CHKERRQ(ierr); 1127 ierr = PetscSectionGetOffset(pointCounts, v, &ooff);CHKERRQ(ierr); 1128 if (dof <= 0) continue; 1129 1130 for (k = 0; k < patch->nsubspaces; ++k) { 1131 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1132 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1133 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1134 PetscInt bs = patch->bs[k]; 1135 PetscInt goff; 1136 1137 for (i = off; i < off + dof; ++i) { 1138 /* Reconstruct mapping of global-to-local on this patch. */ 1139 const PetscInt c = cellsArray[i]; 1140 PetscInt cell = c; 1141 1142 if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);} 1143 for (j = 0; j < nodesPerCell; ++j) { 1144 for (l = 0; l < bs; ++l) { 1145 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset; 1146 const PetscInt localDof = dofsArray[key]; 1147 if (localDof >= 0) {ierr = PetscHMapISet(ht, globalDof, localDof);CHKERRQ(ierr);} 1148 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1149 const PetscInt localDofWithArtificial = dofsArrayWithArtificial[key]; 1150 if (localDofWithArtificial >= 0) { 1151 ierr = PetscHMapISet(htWithArtificial, globalDof, localDofWithArtificial);CHKERRQ(ierr); 1152 } 1153 } 1154 key++; 1155 } 1156 } 1157 } 1158 1159 /* Shove it in the output data structure. */ 1160 ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr); 1161 PetscHashIterBegin(ht, hi); 1162 while (!PetscHashIterAtEnd(ht, hi)) { 1163 PetscInt globalDof, localDof; 1164 1165 PetscHashIterGetKey(ht, hi, globalDof); 1166 PetscHashIterGetVal(ht, hi, localDof); 1167 if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof; 1168 PetscHashIterNext(ht, hi); 1169 } 1170 1171 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1172 ierr = PetscSectionGetOffset(gtolCountsWithArtificial, v, &goff);CHKERRQ(ierr); 1173 PetscHashIterBegin(htWithArtificial, hi); 1174 while (!PetscHashIterAtEnd(htWithArtificial, hi)) { 1175 PetscInt globalDof, localDof; 1176 PetscHashIterGetKey(htWithArtificial, hi, globalDof); 1177 PetscHashIterGetVal(htWithArtificial, hi, localDof); 1178 if (globalDof >= 0) globalDofsArrayWithArtificial[goff + localDof] = globalDof; 1179 PetscHashIterNext(htWithArtificial, hi); 1180 } 1181 } 1182 1183 for (p = 0; p < Np; ++p) { 1184 const PetscInt point = pointsArray[ooff + p]; 1185 PetscInt globalDof, localDof; 1186 1187 ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, point, NULL, &globalDof);CHKERRQ(ierr); 1188 ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr); 1189 offsArray[(ooff + p)*Nf + k] = localDof; 1190 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1191 ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr); 1192 offsArrayWithArtificial[(ooff + p)*Nf + k] = localDof; 1193 } 1194 } 1195 } 1196 1197 ierr = PetscHSetIDestroy(&globalBcs);CHKERRQ(ierr); 1198 ierr = PetscHSetIDestroy(&ownedpts);CHKERRQ(ierr); 1199 ierr = PetscHSetIDestroy(&seenpts);CHKERRQ(ierr); 1200 ierr = PetscHSetIDestroy(&owneddofs);CHKERRQ(ierr); 1201 ierr = PetscHSetIDestroy(&seendofs);CHKERRQ(ierr); 1202 ierr = PetscHSetIDestroy(&artificialbcs);CHKERRQ(ierr); 1203 1204 /* At this point, we have a hash table ht built that maps globalDof -> localDof. 1205 We need to create the dof table laid out cellwise first, then by subspace, 1206 as the assembler assembles cell-wise and we need to stuff the different 1207 contributions of the different function spaces to the right places. So we loop 1208 over cells, then over subspaces. */ 1209 if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */ 1210 for (i = off; i < off + dof; ++i) { 1211 const PetscInt c = cellsArray[i]; 1212 PetscInt cell = c; 1213 1214 if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);} 1215 for (k = 0; k < patch->nsubspaces; ++k) { 1216 const PetscInt *cellNodeMap = patch->cellNodeMap[k]; 1217 PetscInt nodesPerCell = patch->nodesPerCell[k]; 1218 PetscInt subspaceOffset = patch->subspaceOffsets[k]; 1219 PetscInt bs = patch->bs[k]; 1220 1221 for (j = 0; j < nodesPerCell; ++j) { 1222 for (l = 0; l < bs; ++l) { 1223 const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset; 1224 PetscInt localDof; 1225 1226 ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr); 1227 /* If it's not in the hash table, i.e. is a BC dof, 1228 then the PetscHSetIMap above gives -1, which matches 1229 exactly the convention for PETSc's matrix assembly to 1230 ignore the dof. So we don't need to do anything here */ 1231 asmArray[asmKey] = localDof; 1232 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1233 ierr = PetscHMapIGet(htWithArtificial, globalDof, &localDof);CHKERRQ(ierr); 1234 asmArrayWithArtificial[asmKey] = localDof; 1235 } 1236 asmKey++; 1237 } 1238 } 1239 } 1240 } 1241 } 1242 } 1243 if (1 == patch->nsubspaces) { 1244 ierr = PetscMemcpy(asmArray, dofsArray, numDofs * sizeof(PetscInt));CHKERRQ(ierr); 1245 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1246 ierr = PetscMemcpy(asmArrayWithArtificial, dofsArrayWithArtificial, numDofs * sizeof(PetscInt));CHKERRQ(ierr); 1247 } 1248 } 1249 1250 ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr); 1251 ierr = PetscHMapIDestroy(&htWithArtificial);CHKERRQ(ierr); 1252 ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr); 1253 ierr = ISRestoreIndices(points, &pointsArray);CHKERRQ(ierr); 1254 ierr = PetscFree(dofsArray);CHKERRQ(ierr); 1255 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1256 ierr = PetscFree(dofsArrayWithArtificial);CHKERRQ(ierr); 1257 } 1258 /* Create placeholder section for map from points to patch dofs */ 1259 ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->patchSection);CHKERRQ(ierr); 1260 ierr = PetscSectionSetNumFields(patch->patchSection, patch->nsubspaces);CHKERRQ(ierr); 1261 if (patch->combined) { 1262 PetscInt numFields; 1263 ierr = PetscSectionGetNumFields(patch->dofSection[0], &numFields);CHKERRQ(ierr); 1264 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); 1265 ierr = PetscSectionGetChart(patch->dofSection[0], &pStart, &pEnd);CHKERRQ(ierr); 1266 ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr); 1267 for (p = pStart; p < pEnd; ++p) { 1268 PetscInt dof, fdof, f; 1269 1270 ierr = PetscSectionGetDof(patch->dofSection[0], p, &dof);CHKERRQ(ierr); 1271 ierr = PetscSectionSetDof(patch->patchSection, p, dof);CHKERRQ(ierr); 1272 for (f = 0; f < patch->nsubspaces; ++f) { 1273 ierr = PetscSectionGetFieldDof(patch->dofSection[0], p, f, &fdof);CHKERRQ(ierr); 1274 ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr); 1275 } 1276 } 1277 } else { 1278 PetscInt pStartf, pEndf, f; 1279 pStart = PETSC_MAX_INT; 1280 pEnd = PETSC_MIN_INT; 1281 for (f = 0; f < patch->nsubspaces; ++f) { 1282 ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr); 1283 pStart = PetscMin(pStart, pStartf); 1284 pEnd = PetscMax(pEnd, pEndf); 1285 } 1286 ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr); 1287 for (f = 0; f < patch->nsubspaces; ++f) { 1288 ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr); 1289 for (p = pStartf; p < pEndf; ++p) { 1290 PetscInt fdof; 1291 ierr = PetscSectionGetDof(patch->dofSection[f], p, &fdof);CHKERRQ(ierr); 1292 ierr = PetscSectionAddDof(patch->patchSection, p, fdof);CHKERRQ(ierr); 1293 ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr); 1294 } 1295 } 1296 } 1297 ierr = PetscSectionSetUp(patch->patchSection);CHKERRQ(ierr); 1298 ierr = PetscSectionSetUseFieldOffsets(patch->patchSection, PETSC_TRUE);CHKERRQ(ierr); 1299 /* Replace cell indices with firedrake-numbered ones. */ 1300 ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr); 1301 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr); 1302 ierr = PetscObjectSetName((PetscObject) patch->gtol, "Global Indices");CHKERRQ(ierr); 1303 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_g2l_view", patch->classname);CHKERRQ(ierr); 1304 ierr = PetscSectionViewFromOptions(patch->gtolCounts, (PetscObject) pc, option);CHKERRQ(ierr); 1305 ierr = ISViewFromOptions(patch->gtol, (PetscObject) pc, option);CHKERRQ(ierr); 1306 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr); 1307 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArray, PETSC_OWN_POINTER, &patch->offs);CHKERRQ(ierr); 1308 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1309 ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofsWithArtificial, globalDofsArrayWithArtificial, PETSC_OWN_POINTER, &patch->gtolWithArtificial);CHKERRQ(ierr); 1310 ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArrayWithArtificial, PETSC_OWN_POINTER, &patch->dofsWithArtificial);CHKERRQ(ierr); 1311 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArrayWithArtificial, PETSC_OWN_POINTER, &patch->offsWithArtificial);CHKERRQ(ierr); 1312 } 1313 PetscFunctionReturn(0); 1314 } 1315 1316 static PetscErrorCode PCPatchZeroFillMatrix_Private(Mat mat, const PetscInt ncell, const PetscInt ndof, const PetscInt *dof) 1317 { 1318 PetscScalar *values = NULL; 1319 PetscInt rows, c, i; 1320 PetscErrorCode ierr; 1321 1322 PetscFunctionBegin; 1323 ierr = PetscCalloc1(ndof*ndof, &values);CHKERRQ(ierr); 1324 for (c = 0; c < ncell; ++c) { 1325 const PetscInt *idx = &dof[ndof*c]; 1326 ierr = MatSetValues(mat, ndof, idx, ndof, idx, values, INSERT_VALUES);CHKERRQ(ierr); 1327 } 1328 ierr = MatGetLocalSize(mat, &rows, NULL);CHKERRQ(ierr); 1329 for (i = 0; i < rows; ++i) { 1330 ierr = MatSetValues(mat, 1, &i, 1, &i, values, INSERT_VALUES);CHKERRQ(ierr); 1331 } 1332 ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1333 ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1334 ierr = PetscFree(values);CHKERRQ(ierr); 1335 PetscFunctionReturn(0); 1336 } 1337 1338 static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat, PetscBool withArtificial) 1339 { 1340 PC_PATCH *patch = (PC_PATCH *) pc->data; 1341 Vec x, y; 1342 PetscBool flg; 1343 PetscInt csize, rsize; 1344 const char *prefix = NULL; 1345 PetscErrorCode ierr; 1346 1347 PetscFunctionBegin; 1348 if(withArtificial) { 1349 /* would be nice if we could create a rectangular matrix of size numDofsWithArtificial x numDofs here */ 1350 x = patch->patchRHSWithArtificial[point]; 1351 y = patch->patchRHSWithArtificial[point]; 1352 } else { 1353 x = patch->patchRHS[point]; 1354 y = patch->patchUpdate[point]; 1355 } 1356 1357 ierr = VecGetSize(x, &csize);CHKERRQ(ierr); 1358 ierr = VecGetSize(y, &rsize);CHKERRQ(ierr); 1359 ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr); 1360 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 1361 ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr); 1362 ierr = MatAppendOptionsPrefix(*mat, "pc_patch_sub_");CHKERRQ(ierr); 1363 if (patch->sub_mat_type) {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);} 1364 else if (!patch->sub_mat_type) {ierr = MatSetType(*mat, MATDENSE);CHKERRQ(ierr);} 1365 ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr); 1366 ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr); 1367 if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);} 1368 /* Sparse patch matrices */ 1369 if (!flg) { 1370 PetscBT bt; 1371 PetscInt *dnnz = NULL; 1372 const PetscInt *dofsArray = NULL; 1373 PetscInt pStart, pEnd, ncell, offset, c, i, j; 1374 1375 if(withArtificial) { 1376 ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 1377 } else { 1378 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1379 } 1380 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1381 point += pStart; 1382 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1383 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1384 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1385 ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr); 1386 ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1387 /* XXX: This uses N^2 bits to store the sparsity pattern on a 1388 * patch. This is probably OK if the patches are not too big, 1389 * but could use quite a bit of memory for planes in 3D. 1390 * Should we switch based on the value of rsize to a 1391 * hash-table (slower, but more memory efficient) approach? */ 1392 ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr); 1393 for (c = 0; c < ncell; ++c) { 1394 const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell; 1395 for (i = 0; i < patch->totalDofsPerCell; ++i) { 1396 const PetscInt row = idx[i]; 1397 if (row < 0) continue; 1398 for (j = 0; j < patch->totalDofsPerCell; ++j) { 1399 const PetscInt col = idx[j]; 1400 const PetscInt key = row*rsize + col; 1401 if (col < 0) continue; 1402 if (!PetscBTLookupSet(bt, key)) ++dnnz[row]; 1403 } 1404 } 1405 } 1406 ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 1407 ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr); 1408 ierr = PetscFree(dnnz);CHKERRQ(ierr); 1409 ierr = PCPatchZeroFillMatrix_Private(*mat, ncell, patch->totalDofsPerCell, &dofsArray[offset*patch->totalDofsPerCell]);CHKERRQ(ierr); 1410 ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr); 1411 if(withArtificial) { 1412 ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 1413 } else { 1414 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1415 } 1416 } 1417 ierr = MatSetUp(*mat);CHKERRQ(ierr); 1418 PetscFunctionReturn(0); 1419 } 1420 1421 static PetscErrorCode PCPatchComputeFunction_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Vec F, IS cellIS, PetscInt n, const PetscInt *l2p, const PetscInt *l2pWithArtificial, void *ctx) 1422 { 1423 PC_PATCH *patch = (PC_PATCH *) pc->data; 1424 DM dm; 1425 PetscSection s; 1426 const PetscInt *parray, *oarray; 1427 PetscInt Nf = patch->nsubspaces, Np, poff, p, f; 1428 PetscErrorCode ierr; 1429 1430 PetscFunctionBegin; 1431 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 1432 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 1433 /* Set offset into patch */ 1434 ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr); 1435 ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr); 1436 ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr); 1437 ierr = ISGetIndices(patch->offs, &oarray);CHKERRQ(ierr); 1438 for (f = 0; f < Nf; ++f) { 1439 for (p = 0; p < Np; ++p) { 1440 const PetscInt point = parray[poff+p]; 1441 PetscInt dof; 1442 1443 ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr); 1444 ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr); 1445 if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);} 1446 else {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);} 1447 } 1448 } 1449 ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr); 1450 ierr = ISRestoreIndices(patch->offs, &oarray);CHKERRQ(ierr); 1451 if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);} 1452 ierr = DMPlexComputeResidual_Patch_Internal(pc->dm, patch->patchSection, cellIS, 0.0, x, NULL, F, ctx);CHKERRQ(ierr); 1453 PetscFunctionReturn(0); 1454 } 1455 1456 PetscErrorCode PCPatchComputeFunction_Internal(PC pc, Vec x, Vec F, PetscInt point) 1457 { 1458 PC_PATCH *patch = (PC_PATCH *) pc->data; 1459 const PetscInt *dofsArray; 1460 const PetscInt *dofsArrayWithArtificial; 1461 const PetscInt *cellsArray; 1462 PetscInt ncell, offset, pStart, pEnd; 1463 PetscErrorCode ierr; 1464 1465 PetscFunctionBegin; 1466 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1467 if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n"); 1468 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1469 ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArrayWithArtificial);CHKERRQ(ierr); 1470 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1471 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1472 1473 point += pStart; 1474 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1475 1476 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1477 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1478 if (ncell <= 0) { 1479 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1480 PetscFunctionReturn(0); 1481 } 1482 PetscStackPush("PCPatch user callback"); 1483 /* Cannot reuse the same IS because the geometry info is being cached in it */ 1484 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr); 1485 ierr = patch->usercomputef(pc, point, x, F, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, 1486 dofsArrayWithArtificial + offset*patch->totalDofsPerCell, 1487 patch->usercomputefctx);CHKERRQ(ierr); 1488 PetscStackPop; 1489 ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr); 1490 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1491 ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArrayWithArtificial);CHKERRQ(ierr); 1492 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1493 if (patch->viewMatrix) { 1494 char name[PETSC_MAX_PATH_LEN]; 1495 1496 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch vector for Point %D", point);CHKERRQ(ierr); 1497 ierr = PetscObjectSetName((PetscObject) F, name);CHKERRQ(ierr); 1498 ierr = ObjectView((PetscObject) F, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr); 1499 } 1500 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1501 PetscFunctionReturn(0); 1502 } 1503 1504 static PetscErrorCode PCPatchComputeOperator_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Mat J, IS cellIS, PetscInt n, const PetscInt *l2p, const PetscInt *l2pWithArtificial, void *ctx) 1505 { 1506 PC_PATCH *patch = (PC_PATCH *) pc->data; 1507 DM dm; 1508 PetscSection s; 1509 const PetscInt *parray, *oarray; 1510 PetscInt Nf = patch->nsubspaces, Np, poff, p, f; 1511 PetscErrorCode ierr; 1512 1513 PetscFunctionBegin; 1514 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 1515 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 1516 /* Set offset into patch */ 1517 ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr); 1518 ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr); 1519 ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr); 1520 ierr = ISGetIndices(patch->offs, &oarray);CHKERRQ(ierr); 1521 for (f = 0; f < Nf; ++f) { 1522 for (p = 0; p < Np; ++p) { 1523 const PetscInt point = parray[poff+p]; 1524 PetscInt dof; 1525 1526 ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr); 1527 ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr); 1528 if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);} 1529 else {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);} 1530 } 1531 } 1532 ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr); 1533 ierr = ISRestoreIndices(patch->offs, &oarray);CHKERRQ(ierr); 1534 if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);} 1535 /* TODO Shut off MatViewFromOptions() in MatAssemblyEnd() here */ 1536 ierr = DMPlexComputeJacobian_Patch_Internal(pc->dm, patch->patchSection, patch->patchSection, cellIS, 0.0, 0.0, x, NULL, J, J, ctx);CHKERRQ(ierr); 1537 PetscFunctionReturn(0); 1538 } 1539 1540 PetscErrorCode PCPatchComputeOperator_Internal(PC pc, Vec x, Mat mat, PetscInt point, PetscBool withArtificial) 1541 { 1542 PC_PATCH *patch = (PC_PATCH *) pc->data; 1543 const PetscInt *dofsArray; 1544 const PetscInt *dofsArrayWithArtificial = NULL; 1545 const PetscInt *cellsArray; 1546 PetscInt ncell, offset, pStart, pEnd; 1547 PetscBool isNonlinear; 1548 PetscErrorCode ierr; 1549 1550 PetscFunctionBegin; 1551 ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1552 ierr = PetscStrncmp(patch->classname, "snes", sizeof("snes"), &isNonlinear); 1553 if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n"); 1554 if(withArtificial) { 1555 ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 1556 } else { 1557 ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1558 } 1559 if (isNonlinear) { 1560 ierr = ISGetIndices(patch->dofsWithArtificial, &dofsArrayWithArtificial);CHKERRQ(ierr); 1561 } 1562 ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1563 ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr); 1564 1565 point += pStart; 1566 if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr); 1567 1568 ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr); 1569 ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr); 1570 if (ncell <= 0) { 1571 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1572 PetscFunctionReturn(0); 1573 } 1574 PetscStackPush("PCPatch user callback"); 1575 /* Cannot reuse the same IS because the geometry info is being cached in it */ 1576 ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr); 1577 ierr = patch->usercomputeop(pc, point, x, mat, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, dofsArrayWithArtificial ? dofsArrayWithArtificial + offset*patch->totalDofsPerCell : NULL, patch->usercomputeopctx);CHKERRQ(ierr); 1578 PetscStackPop; 1579 ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr); 1580 if(withArtificial) { 1581 ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArray);CHKERRQ(ierr); 1582 } else { 1583 ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr); 1584 } 1585 if (isNonlinear) { 1586 ierr = ISRestoreIndices(patch->dofsWithArtificial, &dofsArrayWithArtificial);CHKERRQ(ierr); 1587 } 1588 ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr); 1589 if (patch->viewMatrix) { 1590 char name[PETSC_MAX_PATH_LEN]; 1591 1592 ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch matrix for Point %D", point);CHKERRQ(ierr); 1593 ierr = PetscObjectSetName((PetscObject) mat, name);CHKERRQ(ierr); 1594 ierr = ObjectView((PetscObject) mat, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr); 1595 } 1596 ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr); 1597 PetscFunctionReturn(0); 1598 } 1599 1600 PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat, PetscBool withArtificial) 1601 { 1602 PC_PATCH *patch = (PC_PATCH *) pc->data; 1603 const PetscScalar *xArray = NULL; 1604 PetscScalar *yArray = NULL; 1605 const PetscInt *gtolArray = NULL; 1606 PetscInt dof, offset, lidx; 1607 PetscErrorCode ierr; 1608 1609 PetscFunctionBeginHot; 1610 ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 1611 ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr); 1612 ierr = VecGetArray(y, &yArray);CHKERRQ(ierr); 1613 if(withArtificial) { 1614 ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr); 1615 ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offset);CHKERRQ(ierr); 1616 ierr = ISGetIndices(patch->gtolWithArtificial, >olArray);CHKERRQ(ierr); 1617 } else { 1618 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 1619 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 1620 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 1621 } 1622 if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n"); 1623 if (mode == ADD_VALUES && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n"); 1624 for (lidx = 0; lidx < dof; ++lidx) { 1625 const PetscInt gidx = gtolArray[offset+lidx]; 1626 1627 if (mode == INSERT_VALUES) yArray[lidx] = xArray[gidx]; /* Forward */ 1628 else yArray[gidx] += xArray[lidx]; /* Reverse */ 1629 } 1630 if(withArtificial) { 1631 ierr = ISRestoreIndices(patch->gtolWithArtificial, >olArray);CHKERRQ(ierr); 1632 } else { 1633 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 1634 } 1635 ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr); 1636 ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr); 1637 ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr); 1638 PetscFunctionReturn(0); 1639 } 1640 1641 static PetscErrorCode PCSetUp_PATCH_Linear(PC pc) 1642 { 1643 PC_PATCH *patch = (PC_PATCH *) pc->data; 1644 const char *prefix; 1645 PetscInt i; 1646 PetscErrorCode ierr; 1647 1648 PetscFunctionBegin; 1649 if (!pc->setupcalled) { 1650 ierr = PetscMalloc1(patch->npatch, &patch->solver);CHKERRQ(ierr); 1651 ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 1652 for (i = 0; i < patch->npatch; ++i) { 1653 KSP ksp; 1654 PC subpc; 1655 1656 ierr = KSPCreate(PETSC_COMM_SELF, &ksp);CHKERRQ(ierr); 1657 ierr = KSPSetOptionsPrefix(ksp, prefix);CHKERRQ(ierr); 1658 ierr = KSPAppendOptionsPrefix(ksp, "sub_");CHKERRQ(ierr); 1659 ierr = PetscObjectIncrementTabLevel((PetscObject) ksp, (PetscObject) pc, 1);CHKERRQ(ierr); 1660 ierr = KSPGetPC(ksp, &subpc);CHKERRQ(ierr); 1661 ierr = PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1);CHKERRQ(ierr); 1662 ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) ksp);CHKERRQ(ierr); 1663 patch->solver[i] = (PetscObject) ksp; 1664 } 1665 } 1666 if (patch->save_operators) { 1667 for (i = 0; i < patch->npatch; ++i) { 1668 ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr); 1669 ierr = PCPatchComputeOperator_Internal(pc, NULL, patch->mat[i], i, PETSC_FALSE);CHKERRQ(ierr); 1670 ierr = KSPSetOperators((KSP) patch->solver[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr); 1671 } 1672 } 1673 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1674 for (i = 0; i < patch->npatch; ++i) { 1675 /* Instead of padding patch->patchUpdate with zeros to get */ 1676 /* patch->patchUpdateWithArtificial and then multiplying with the matrix, */ 1677 /* just get rid of the columns that correspond to the dofs with */ 1678 /* artificial bcs. That's of course fairly inefficient, hopefully we */ 1679 /* can just assemble the rectangular matrix in the first place. */ 1680 Mat matSquare; 1681 IS rowis; 1682 PetscInt dof; 1683 1684 ierr = MatGetSize(patch->mat[i], &dof, NULL);CHKERRQ(ierr); 1685 if (dof == 0) { 1686 patch->matWithArtificial[i] = NULL; 1687 continue; 1688 } 1689 1690 ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr); 1691 ierr = MatZeroEntries(matSquare);CHKERRQ(ierr); 1692 ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr); 1693 1694 ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr); 1695 ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr); 1696 if(pc->setupcalled) { 1697 ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_REUSE_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr); 1698 } else { 1699 ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &patch->matWithArtificial[i]); CHKERRQ(ierr); 1700 } 1701 ierr = ISDestroy(&rowis); CHKERRQ(ierr); 1702 ierr = MatDestroy(&matSquare);CHKERRQ(ierr); 1703 } 1704 } 1705 PetscFunctionReturn(0); 1706 } 1707 1708 static PetscErrorCode PCSetUp_PATCH(PC pc) 1709 { 1710 PC_PATCH *patch = (PC_PATCH *) pc->data; 1711 PetscInt i; 1712 PetscBool isNonlinear; 1713 PetscErrorCode ierr; 1714 1715 PetscFunctionBegin; 1716 if (!pc->setupcalled) { 1717 PetscInt pStart, pEnd, p; 1718 PetscInt localSize; 1719 1720 ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 1721 1722 ierr = PetscStrncmp(patch->classname, "snes", sizeof("snes"), &isNonlinear); 1723 if (!patch->nsubspaces) { 1724 DM dm; 1725 PetscDS prob; 1726 PetscSection s; 1727 PetscInt cStart, cEnd, c, Nf, f, numGlobalBcs = 0, *globalBcs, *Nb, totNb = 0, **cellDofs; 1728 1729 ierr = PCGetDM(pc, &dm);CHKERRQ(ierr); 1730 if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Must set DM for PCPATCH or call PCPatchSetDiscretisationInfo()"); 1731 ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr); 1732 ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr); 1733 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 1734 for (p = pStart; p < pEnd; ++p) { 1735 PetscInt cdof; 1736 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1737 numGlobalBcs += cdof; 1738 } 1739 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1740 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 1741 ierr = PetscMalloc3(Nf, &Nb, Nf, &cellDofs, numGlobalBcs, &globalBcs);CHKERRQ(ierr); 1742 for (f = 0; f < Nf; ++f) { 1743 PetscFE fe; 1744 PetscDualSpace sp; 1745 PetscInt cdoff = 0; 1746 1747 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1748 /* ierr = PetscFEGetNumComponents(fe, &Nc[f]);CHKERRQ(ierr); */ 1749 ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 1750 ierr = PetscDualSpaceGetDimension(sp, &Nb[f]);CHKERRQ(ierr); 1751 totNb += Nb[f]; 1752 1753 ierr = PetscMalloc1((cEnd-cStart)*Nb[f], &cellDofs[f]);CHKERRQ(ierr); 1754 for (c = cStart; c < cEnd; ++c) { 1755 PetscInt *closure = NULL; 1756 PetscInt clSize = 0, cl; 1757 1758 ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 1759 for (cl = 0; cl < clSize*2; cl += 2) { 1760 const PetscInt p = closure[cl]; 1761 PetscInt fdof, d, foff; 1762 1763 ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr); 1764 ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr); 1765 for (d = 0; d < fdof; ++d, ++cdoff) cellDofs[f][cdoff] = foff + d; 1766 } 1767 ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 1768 } 1769 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]); 1770 } 1771 numGlobalBcs = 0; 1772 for (p = pStart; p < pEnd; ++p) { 1773 const PetscInt *ind; 1774 PetscInt off, cdof, d; 1775 1776 ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr); 1777 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 1778 ierr = PetscSectionGetConstraintIndices(s, p, &ind);CHKERRQ(ierr); 1779 for (d = 0; d < cdof; ++d) globalBcs[numGlobalBcs++] = off + ind[d]; 1780 } 1781 1782 ierr = PCPatchSetDiscretisationInfoCombined(pc, dm, Nb, (const PetscInt **) cellDofs, numGlobalBcs, globalBcs, numGlobalBcs, globalBcs);CHKERRQ(ierr); 1783 for (f = 0; f < Nf; ++f) { 1784 ierr = PetscFree(cellDofs[f]);CHKERRQ(ierr); 1785 } 1786 ierr = PetscFree3(Nb, cellDofs, globalBcs);CHKERRQ(ierr); 1787 ierr = PCPatchSetComputeFunction(pc, PCPatchComputeFunction_DMPlex_Private, NULL);CHKERRQ(ierr); 1788 ierr = PCPatchSetComputeOperator(pc, PCPatchComputeOperator_DMPlex_Private, NULL);CHKERRQ(ierr); 1789 } 1790 1791 localSize = patch->subspaceOffsets[patch->nsubspaces]; 1792 ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localRHS);CHKERRQ(ierr); 1793 ierr = VecSetUp(patch->localRHS);CHKERRQ(ierr); 1794 ierr = VecDuplicate(patch->localRHS, &patch->localUpdate);CHKERRQ(ierr); 1795 ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr); 1796 ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr); 1797 1798 /* OK, now build the work vectors */ 1799 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr); 1800 ierr = PetscMalloc1(patch->npatch, &patch->patchRHS);CHKERRQ(ierr); 1801 ierr = PetscMalloc1(patch->npatch, &patch->patchUpdate);CHKERRQ(ierr); 1802 1803 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 1804 ierr = PetscMalloc1(patch->npatch, &patch->patchRHSWithArtificial);CHKERRQ(ierr); 1805 } 1806 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1807 ierr = PetscMalloc1(patch->npatch, &patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr); 1808 } 1809 for (p = pStart; p < pEnd; ++p) { 1810 PetscInt dof; 1811 1812 ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr); 1813 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHS[p-pStart]);CHKERRQ(ierr); 1814 ierr = VecSetUp(patch->patchRHS[p-pStart]);CHKERRQ(ierr); 1815 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchUpdate[p-pStart]);CHKERRQ(ierr); 1816 ierr = VecSetUp(patch->patchUpdate[p-pStart]);CHKERRQ(ierr); 1817 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE || isNonlinear) { 1818 const PetscInt *gtolArray, *gtolArrayWithArtificial = NULL; 1819 PetscInt numPatchDofs, offset; 1820 PetscInt numPatchDofsWithArtificial, offsetWithArtificial; 1821 PetscInt dofWithoutArtificialCounter = 0; 1822 PetscInt *patchWithoutArtificialToWithArtificialArray; 1823 1824 if (patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { /* no || isNonlinear */ 1825 ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &dof);CHKERRQ(ierr); 1826 ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr); 1827 ierr = VecSetUp(patch->patchRHSWithArtificial[p-pStart]);CHKERRQ(ierr); 1828 } 1829 1830 /* Now build the mapping that for a dof in a patch WITHOUT dofs that have artificial bcs gives the */ 1831 /* the index in the patch with all dofs */ 1832 ierr = ISGetIndices(patch->gtol, >olArray);CHKERRQ(ierr); 1833 1834 ierr = PetscSectionGetDof(patch->gtolCounts, p, &numPatchDofs);CHKERRQ(ierr); 1835 1836 ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr); 1837 ierr = ISGetIndices(patch->gtolWithArtificial, >olArrayWithArtificial);CHKERRQ(ierr); 1838 ierr = PetscSectionGetDof(patch->gtolCountsWithArtificial, p, &numPatchDofsWithArtificial);CHKERRQ(ierr); 1839 ierr = PetscSectionGetOffset(patch->gtolCountsWithArtificial, p, &offsetWithArtificial);CHKERRQ(ierr); 1840 1841 ierr = PetscMalloc1(numPatchDofs, &patchWithoutArtificialToWithArtificialArray);CHKERRQ(ierr); 1842 1843 ierr = ISCreateGeneral(PETSC_COMM_SELF, numPatchDofs, patchWithoutArtificialToWithArtificialArray, PETSC_OWN_POINTER, &patch->dofMappingWithoutToWithArtificial[p-pStart]);CHKERRQ(ierr); 1844 if (numPatchDofs == 0) continue; 1845 for (i=0; i<numPatchDofsWithArtificial; i++) { 1846 if (gtolArrayWithArtificial[i+offsetWithArtificial] == gtolArray[offset+dofWithoutArtificialCounter]) { 1847 patchWithoutArtificialToWithArtificialArray[dofWithoutArtificialCounter] = i; 1848 dofWithoutArtificialCounter++; 1849 if (dofWithoutArtificialCounter == numPatchDofs) 1850 break; 1851 } 1852 } 1853 ierr = ISRestoreIndices(patch->gtol, >olArray);CHKERRQ(ierr); 1854 ierr = ISRestoreIndices(patch->gtolWithArtificial, >olArrayWithArtificial);CHKERRQ(ierr); 1855 } 1856 } 1857 if (patch->save_operators) { 1858 ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr); 1859 for (i = 0; i < patch->npatch; ++i) { 1860 ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i], PETSC_FALSE);CHKERRQ(ierr); 1861 } 1862 } 1863 ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr); 1864 1865 /* If desired, calculate weights for dof multiplicity */ 1866 if (patch->partition_of_unity) { 1867 PetscScalar *input = NULL; 1868 PetscScalar *output = NULL; 1869 Vec global; 1870 1871 ierr = VecDuplicate(patch->localRHS, &patch->dof_weights);CHKERRQ(ierr); 1872 if(patch->local_composition_type == PC_COMPOSITE_ADDITIVE) { 1873 for (i = 0; i < patch->npatch; ++i) { 1874 PetscInt dof; 1875 1876 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr); 1877 if (dof <= 0) continue; 1878 ierr = VecSet(patch->patchRHS[i], 1.0);CHKERRQ(ierr); 1879 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchRHS[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE, PETSC_FALSE);CHKERRQ(ierr); 1880 } 1881 } else { 1882 /* multiplicative is actually only locally multiplicative and globally additive. need the pou where the mesh decomposition overlaps */ 1883 ierr = VecSet(patch->dof_weights, 1.0);CHKERRQ(ierr); 1884 } 1885 1886 VecDuplicate(patch->dof_weights, &global); 1887 VecSet(global, 0.); 1888 1889 ierr = VecGetArray(patch->dof_weights, &input);CHKERRQ(ierr); 1890 ierr = VecGetArray(global, &output);CHKERRQ(ierr); 1891 ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr); 1892 ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, input, output, MPI_SUM);CHKERRQ(ierr); 1893 ierr = VecRestoreArray(patch->dof_weights, &input);CHKERRQ(ierr); 1894 ierr = VecRestoreArray(global, &output);CHKERRQ(ierr); 1895 1896 ierr = VecReciprocal(global);CHKERRQ(ierr); 1897 1898 ierr = VecGetArray(patch->dof_weights, &output);CHKERRQ(ierr); 1899 ierr = VecGetArray(global, &input);CHKERRQ(ierr); 1900 ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr); 1901 ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, input, output);CHKERRQ(ierr); 1902 ierr = VecRestoreArray(patch->dof_weights, &output);CHKERRQ(ierr); 1903 ierr = VecRestoreArray(global, &input);CHKERRQ(ierr); 1904 ierr = VecDestroy(&global);CHKERRQ(ierr); 1905 } 1906 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE && patch->save_operators) { 1907 ierr = PetscMalloc1(patch->npatch, &patch->matWithArtificial);CHKERRQ(ierr); 1908 } 1909 } 1910 ierr = (*patch->setupsolver)(pc);CHKERRQ(ierr); 1911 PetscFunctionReturn(0); 1912 } 1913 1914 static PetscErrorCode PCApply_PATCH_Linear(PC pc, PetscInt i, Vec x, Vec y) 1915 { 1916 PC_PATCH *patch = (PC_PATCH *) pc->data; 1917 KSP ksp = (KSP) patch->solver[i]; 1918 PetscErrorCode ierr; 1919 1920 PetscFunctionBegin; 1921 if (!patch->save_operators) { 1922 Mat mat; 1923 1924 ierr = PCPatchCreateMatrix_Private(pc, i, &mat, PETSC_FALSE);CHKERRQ(ierr); 1925 /* Populate operator here. */ 1926 ierr = PCPatchComputeOperator_Internal(pc, NULL, mat, i, PETSC_FALSE);CHKERRQ(ierr); 1927 ierr = KSPSetOperators(ksp, mat, mat);CHKERRQ(ierr); 1928 /* Drop reference so the KSPSetOperators below will blow it away. */ 1929 ierr = MatDestroy(&mat);CHKERRQ(ierr); 1930 } 1931 ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 1932 if (!ksp->setfromoptionscalled) { 1933 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 1934 } 1935 ierr = KSPSolve(ksp, x, y);CHKERRQ(ierr); 1936 ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr); 1937 if (!patch->save_operators) { 1938 PC pc; 1939 ierr = KSPSetOperators(ksp, NULL, NULL);CHKERRQ(ierr); 1940 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 1941 /* Destroy PC context too, otherwise the factored matrix hangs around. */ 1942 ierr = PCReset(pc);CHKERRQ(ierr); 1943 } 1944 PetscFunctionReturn(0); 1945 } 1946 1947 static PetscErrorCode PCUpdateMultiplicative_PATCH_Linear(PC pc, PetscInt i, PetscInt pStart) 1948 { 1949 PC_PATCH *patch = (PC_PATCH *) pc->data; 1950 Mat multMat; 1951 PetscErrorCode ierr; 1952 1953 PetscFunctionBegin; 1954 1955 if (patch->save_operators) { 1956 multMat = patch->matWithArtificial[i]; 1957 } else { 1958 /*Very inefficient, hopefully we can just assemble the rectangular matrix in the first place.*/ 1959 Mat matSquare; 1960 PetscInt dof; 1961 IS rowis; 1962 ierr = PCPatchCreateMatrix_Private(pc, i, &matSquare, PETSC_TRUE);CHKERRQ(ierr); 1963 ierr = MatZeroEntries(matSquare);CHKERRQ(ierr); 1964 ierr = PCPatchComputeOperator_Internal(pc, NULL, matSquare, i, PETSC_TRUE);CHKERRQ(ierr); 1965 ierr = MatGetSize(matSquare, &dof, NULL);CHKERRQ(ierr); 1966 ierr = ISCreateStride(PETSC_COMM_SELF, dof, 0, 1, &rowis); CHKERRQ(ierr); 1967 ierr = MatCreateSubMatrix(matSquare, rowis, patch->dofMappingWithoutToWithArtificial[i], MAT_INITIAL_MATRIX, &multMat); CHKERRQ(ierr); 1968 ierr = MatDestroy(&matSquare);CHKERRQ(ierr); 1969 ierr = ISDestroy(&rowis); CHKERRQ(ierr); 1970 } 1971 ierr = MatMult(multMat, patch->patchUpdate[i], patch->patchRHSWithArtificial[i]); CHKERRQ(ierr); 1972 ierr = VecScale(patch->patchRHSWithArtificial[i], -1.0); CHKERRQ(ierr); 1973 ierr = PCPatch_ScatterLocal_Private(pc, i + pStart, patch->patchRHSWithArtificial[i], patch->localRHS, ADD_VALUES, SCATTER_REVERSE, PETSC_TRUE); CHKERRQ(ierr); 1974 if (!patch->save_operators) { 1975 ierr = MatDestroy(&multMat); CHKERRQ(ierr); 1976 } 1977 PetscFunctionReturn(0); 1978 } 1979 1980 static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y) 1981 { 1982 PC_PATCH *patch = (PC_PATCH *) pc->data; 1983 const PetscScalar *globalRHS = NULL; 1984 PetscScalar *localRHS = NULL; 1985 PetscScalar *globalUpdate = NULL; 1986 const PetscInt *bcNodes = NULL; 1987 PetscInt nsweep = patch->symmetrise_sweep ? 2 : 1; 1988 PetscInt start[2] = {0, 0}; 1989 PetscInt end[2] = {-1, -1}; 1990 const PetscInt inc[2] = {1, -1}; 1991 const PetscScalar *localUpdate; 1992 const PetscInt *iterationSet; 1993 PetscInt pStart, numBcs, n, sweep, bc, j; 1994 PetscErrorCode ierr; 1995 1996 PetscFunctionBegin; 1997 ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 1998 ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr); 1999 /* start, end, inc have 2 entries to manage a second backward sweep if we symmetrize */ 2000 end[0] = patch->npatch; 2001 start[1] = patch->npatch-1; 2002 if (patch->user_patches) { 2003 ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr); 2004 start[1] = end[0] - 1; 2005 ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr); 2006 } 2007 /* Scatter from global space into overlapped local spaces */ 2008 ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr); 2009 ierr = VecGetArray(patch->localRHS, &localRHS);CHKERRQ(ierr); 2010 ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr); 2011 ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalRHS, localRHS);CHKERRQ(ierr); 2012 ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr); 2013 ierr = VecRestoreArray(patch->localRHS, &localRHS);CHKERRQ(ierr); 2014 2015 ierr = VecSet(patch->localUpdate, 0.0);CHKERRQ(ierr); 2016 ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr); 2017 for (sweep = 0; sweep < nsweep; sweep++) { 2018 for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) { 2019 PetscInt i = patch->user_patches ? iterationSet[j] : j; 2020 PetscInt start, len; 2021 2022 ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr); 2023 ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr); 2024 /* TODO: Squash out these guys in the setup as well. */ 2025 if (len <= 0) continue; 2026 /* TODO: Do we need different scatters for X and Y? */ 2027 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localRHS, patch->patchRHS[i], INSERT_VALUES, SCATTER_FORWARD, PETSC_FALSE);CHKERRQ(ierr); 2028 ierr = (*patch->applysolver)(pc, i, patch->patchRHS[i], patch->patchUpdate[i]);CHKERRQ(ierr); 2029 ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchUpdate[i], patch->localUpdate, ADD_VALUES, SCATTER_REVERSE, PETSC_FALSE);CHKERRQ(ierr); 2030 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2031 ierr = (*patch->updatemultiplicative)(pc, i, pStart);CHKERRQ(ierr); 2032 } 2033 } 2034 } 2035 if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);} 2036 /* XXX: should we do this on the global vector? */ 2037 if (patch->partition_of_unity) { 2038 ierr = VecPointwiseMult(patch->localUpdate, patch->localUpdate, patch->dof_weights);CHKERRQ(ierr); 2039 } 2040 /* Now patch->localUpdate contains the solution of the patch solves, so we need to combine them all. */ 2041 ierr = VecSet(y, 0.0);CHKERRQ(ierr); 2042 ierr = VecGetArray(y, &globalUpdate);CHKERRQ(ierr); 2043 ierr = VecGetArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr); 2044 ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr); 2045 ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localUpdate, globalUpdate, MPI_SUM);CHKERRQ(ierr); 2046 ierr = VecRestoreArrayRead(patch->localUpdate, &localUpdate);CHKERRQ(ierr); 2047 2048 /* Now we need to send the global BC values through */ 2049 ierr = VecGetArrayRead(x, &globalRHS);CHKERRQ(ierr); 2050 ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr); 2051 ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 2052 ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr); 2053 for (bc = 0; bc < numBcs; ++bc) { 2054 const PetscInt idx = bcNodes[bc]; 2055 if (idx < n) globalUpdate[idx] = globalRHS[idx]; 2056 } 2057 2058 ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr); 2059 ierr = VecRestoreArrayRead(x, &globalRHS);CHKERRQ(ierr); 2060 ierr = VecRestoreArray(y, &globalUpdate);CHKERRQ(ierr); 2061 2062 ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr); 2063 ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr); 2064 PetscFunctionReturn(0); 2065 } 2066 2067 static PetscErrorCode PCReset_PATCH_Linear(PC pc) 2068 { 2069 PC_PATCH *patch = (PC_PATCH *) pc->data; 2070 PetscInt i; 2071 PetscErrorCode ierr; 2072 2073 PetscFunctionBegin; 2074 if (patch->solver) { 2075 for (i = 0; i < patch->npatch; ++i) {ierr = KSPReset((KSP) patch->solver[i]);CHKERRQ(ierr);} 2076 } 2077 PetscFunctionReturn(0); 2078 } 2079 2080 static PetscErrorCode PCReset_PATCH(PC pc) 2081 { 2082 PC_PATCH *patch = (PC_PATCH *) pc->data; 2083 PetscInt i; 2084 PetscErrorCode ierr; 2085 2086 PetscFunctionBegin; 2087 /* TODO: Get rid of all these ifs */ 2088 ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr); 2089 ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr); 2090 ierr = PetscSectionDestroy(&patch->pointCounts);CHKERRQ(ierr); 2091 ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr); 2092 ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr); 2093 ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr); 2094 ierr = ISDestroy(&patch->cells);CHKERRQ(ierr); 2095 ierr = ISDestroy(&patch->points);CHKERRQ(ierr); 2096 ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr); 2097 ierr = ISDestroy(&patch->offs);CHKERRQ(ierr); 2098 ierr = PetscSectionDestroy(&patch->patchSection);CHKERRQ(ierr); 2099 ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr); 2100 ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr); 2101 ierr = PetscSectionDestroy(&patch->gtolCountsWithArtificial);CHKERRQ(ierr); 2102 ierr = ISDestroy(&patch->gtolWithArtificial);CHKERRQ(ierr); 2103 ierr = ISDestroy(&patch->dofsWithArtificial);CHKERRQ(ierr); 2104 ierr = ISDestroy(&patch->offsWithArtificial);CHKERRQ(ierr); 2105 2106 2107 if (patch->dofSection) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);} 2108 ierr = PetscFree(patch->dofSection);CHKERRQ(ierr); 2109 ierr = PetscFree(patch->bs);CHKERRQ(ierr); 2110 ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr); 2111 if (patch->cellNodeMap) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscFree(patch->cellNodeMap[i]);CHKERRQ(ierr);} 2112 ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr); 2113 ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr); 2114 2115 ierr = (*patch->resetsolver)(pc);CHKERRQ(ierr); 2116 2117 if (patch->subspaces_to_exclude) { 2118 PetscHSetIDestroy(&patch->subspaces_to_exclude); 2119 } 2120 2121 ierr = VecDestroy(&patch->localRHS);CHKERRQ(ierr); 2122 ierr = VecDestroy(&patch->localUpdate);CHKERRQ(ierr); 2123 if (patch->patchRHS) { 2124 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHS[i]);CHKERRQ(ierr);} 2125 ierr = PetscFree(patch->patchRHS);CHKERRQ(ierr); 2126 } 2127 if (patch->patchUpdate) { 2128 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchUpdate[i]);CHKERRQ(ierr);} 2129 ierr = PetscFree(patch->patchUpdate);CHKERRQ(ierr); 2130 } 2131 ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr); 2132 if (patch->patch_dof_weights) { 2133 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);} 2134 ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr); 2135 } 2136 if (patch->mat) { 2137 for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);} 2138 ierr = PetscFree(patch->mat);CHKERRQ(ierr); 2139 } 2140 if (patch->matWithArtificial) { 2141 for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->matWithArtificial[i]);CHKERRQ(ierr);} 2142 ierr = PetscFree(patch->matWithArtificial);CHKERRQ(ierr); 2143 } 2144 if (patch->patchRHSWithArtificial) { 2145 for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchRHSWithArtificial[i]);CHKERRQ(ierr);} 2146 ierr = PetscFree(patch->patchRHSWithArtificial);CHKERRQ(ierr); 2147 } 2148 if(patch->dofMappingWithoutToWithArtificial) { 2149 for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->dofMappingWithoutToWithArtificial[i]);CHKERRQ(ierr);} 2150 ierr = PetscFree(patch->dofMappingWithoutToWithArtificial);CHKERRQ(ierr); 2151 2152 } 2153 ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr); 2154 if (patch->userIS) { 2155 for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);} 2156 ierr = PetscFree(patch->userIS);CHKERRQ(ierr); 2157 } 2158 patch->bs = 0; 2159 patch->cellNodeMap = NULL; 2160 patch->nsubspaces = 0; 2161 ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr); 2162 2163 ierr = PetscViewerDestroy(&patch->viewerSection);CHKERRQ(ierr); 2164 PetscFunctionReturn(0); 2165 } 2166 2167 static PetscErrorCode PCDestroy_PATCH_Linear(PC pc) 2168 { 2169 PC_PATCH *patch = (PC_PATCH *) pc->data; 2170 PetscInt i; 2171 PetscErrorCode ierr; 2172 2173 PetscFunctionBegin; 2174 if (patch->solver) { 2175 for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy((KSP *) &patch->solver[i]);CHKERRQ(ierr);} 2176 ierr = PetscFree(patch->solver);CHKERRQ(ierr); 2177 } 2178 PetscFunctionReturn(0); 2179 } 2180 2181 static PetscErrorCode PCDestroy_PATCH(PC pc) 2182 { 2183 PC_PATCH *patch = (PC_PATCH *) pc->data; 2184 PetscErrorCode ierr; 2185 2186 PetscFunctionBegin; 2187 ierr = PCReset_PATCH(pc);CHKERRQ(ierr); 2188 ierr = (*patch->destroysolver)(pc);CHKERRQ(ierr); 2189 ierr = PetscFree(pc->data);CHKERRQ(ierr); 2190 PetscFunctionReturn(0); 2191 } 2192 2193 static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc) 2194 { 2195 PC_PATCH *patch = (PC_PATCH *) pc->data; 2196 PCPatchConstructType patchConstructionType = PC_PATCH_STAR; 2197 char sub_mat_type[PETSC_MAX_PATH_LEN]; 2198 char option[PETSC_MAX_PATH_LEN]; 2199 const char *prefix; 2200 PetscBool flg, dimflg, codimflg; 2201 MPI_Comm comm; 2202 PetscInt *ifields, nfields, k; 2203 PetscErrorCode ierr; 2204 PCCompositeType loctype = PC_COMPOSITE_ADDITIVE; 2205 2206 PetscFunctionBegin; 2207 ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr); 2208 ierr = PetscObjectGetOptionsPrefix((PetscObject) pc, &prefix);CHKERRQ(ierr); 2209 ierr = PetscOptionsHead(PetscOptionsObject, "Patch solver options");CHKERRQ(ierr); 2210 2211 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_save_operators", patch->classname); 2212 ierr = PetscOptionsBool(option, "Store all patch operators for lifetime of object?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr); 2213 2214 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_partition_of_unity", patch->classname); 2215 ierr = PetscOptionsBool(option, "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr); 2216 2217 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_local_type", patch->classname); 2218 ierr = PetscOptionsEnum(option,"Type of local solver composition (additive or multiplicative)","PCPatchSetLocalComposition",PCCompositeTypes,(PetscEnum)loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr); 2219 if(flg) { ierr = PCPatchSetLocalComposition(pc, loctype);CHKERRQ(ierr);} 2220 2221 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_dim", patch->classname); 2222 ierr = PetscOptionsInt(option, "What dimension of mesh point to construct patches by? (0 = vertices)", "PCPATCH", patch->dim, &patch->dim, &dimflg);CHKERRQ(ierr); 2223 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_codim", patch->classname); 2224 ierr = PetscOptionsInt(option, "What co-dimension of mesh point to construct patches by? (0 = cells)", "PCPATCH", patch->codim, &patch->codim, &codimflg);CHKERRQ(ierr); 2225 if (dimflg && codimflg) {SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension");CHKERRQ(ierr);} 2226 2227 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_construct_type", patch->classname); 2228 ierr = PetscOptionsEnum(option, "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr); 2229 if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);} 2230 2231 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_vanka_dim", patch->classname); 2232 ierr = PetscOptionsInt(option, "Topological dimension of entities for Vanka to ignore", "PCPATCH", patch->vankadim, &patch->vankadim, &flg);CHKERRQ(ierr); 2233 2234 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_ignore_dim", patch->classname); 2235 ierr = PetscOptionsInt(option, "Topological dimension of entities for completion to ignore", "PCPATCH", patch->ignoredim, &patch->ignoredim, &flg);CHKERRQ(ierr); 2236 2237 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_sub_mat_type", patch->classname); 2238 ierr = PetscOptionsFList(option, "Matrix type for patch solves", "PCPatchSetSubMatType", MatList, NULL, sub_mat_type, PETSC_MAX_PATH_LEN, &flg);CHKERRQ(ierr); 2239 if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);} 2240 2241 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_symmetrise_sweep", patch->classname); 2242 ierr = PetscOptionsBool(option, "Go start->end, end->start?", "PCPATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr); 2243 2244 /* If the user has set the number of subspaces, use that for the buffer size, 2245 otherwise use a large number */ 2246 if (patch->nsubspaces <= 0) { 2247 nfields = 128; 2248 } else { 2249 nfields = patch->nsubspaces; 2250 } 2251 ierr = PetscMalloc1(nfields, &ifields);CHKERRQ(ierr); 2252 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_exclude_subspaces", patch->classname); 2253 ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,option,ifields,&nfields,&flg);CHKERRQ(ierr); 2254 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"); 2255 if (flg) { 2256 PetscHSetIClear(patch->subspaces_to_exclude); 2257 for (k = 0; k < nfields; k++) { 2258 PetscHSetIAdd(patch->subspaces_to_exclude, ifields[k]); 2259 } 2260 } 2261 ierr = PetscFree(ifields);CHKERRQ(ierr); 2262 2263 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_patches_view", patch->classname); 2264 ierr = PetscOptionsBool(option, "Print out information during patch construction", "PCPATCH", patch->viewPatches, &patch->viewPatches, &flg);CHKERRQ(ierr); 2265 2266 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_cells_view", patch->classname); 2267 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option, &patch->viewerCells, &patch->formatCells, &patch->viewCells);CHKERRQ(ierr); 2268 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_points_view", patch->classname); 2269 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option, &patch->viewerPoints, &patch->formatPoints, &patch->viewPoints);CHKERRQ(ierr); 2270 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_section_view", patch->classname); 2271 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option, &patch->viewerSection, &patch->formatSection, &patch->viewSection);CHKERRQ(ierr); 2272 ierr = PetscSNPrintf(option, PETSC_MAX_PATH_LEN, "-%s_patch_mat_view", patch->classname); 2273 ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, option, &patch->viewerMatrix, &patch->formatMatrix, &patch->viewMatrix);CHKERRQ(ierr); 2274 ierr = PetscOptionsTail();CHKERRQ(ierr); 2275 patch->optionsSet = PETSC_TRUE; 2276 PetscFunctionReturn(0); 2277 } 2278 2279 static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc) 2280 { 2281 PC_PATCH *patch = (PC_PATCH*) pc->data; 2282 KSPConvergedReason reason; 2283 PetscInt i; 2284 PetscErrorCode ierr; 2285 2286 PetscFunctionBegin; 2287 if (!patch->save_operators) { 2288 /* Can't do this here because the sub KSPs don't have an operator attached yet. */ 2289 PetscFunctionReturn(0); 2290 } 2291 for (i = 0; i < patch->npatch; ++i) { 2292 if (!((KSP) patch->solver[i])->setfromoptionscalled) { 2293 ierr = KSPSetFromOptions((KSP) patch->solver[i]);CHKERRQ(ierr); 2294 } 2295 ierr = KSPSetUp((KSP) patch->solver[i]);CHKERRQ(ierr); 2296 ierr = KSPGetConvergedReason((KSP) patch->solver[i], &reason);CHKERRQ(ierr); 2297 if (reason == KSP_DIVERGED_PCSETUP_FAILED) pc->failedreason = PC_SUBPC_ERROR; 2298 } 2299 PetscFunctionReturn(0); 2300 } 2301 2302 static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer) 2303 { 2304 PC_PATCH *patch = (PC_PATCH *) pc->data; 2305 PetscViewer sviewer; 2306 PetscBool isascii; 2307 PetscMPIInt rank; 2308 PetscErrorCode ierr; 2309 2310 PetscFunctionBegin; 2311 /* TODO Redo tabbing with set tbas in new style */ 2312 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 2313 if (!isascii) PetscFunctionReturn(0); 2314 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr); 2315 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 2316 ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr); 2317 if(patch->local_composition_type == PC_COMPOSITE_MULTIPLICATIVE) { 2318 ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: multiplicative\n");CHKERRQ(ierr); 2319 } else { 2320 ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr); 2321 } 2322 if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);} 2323 else {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);} 2324 if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);} 2325 else {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);} 2326 if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);} 2327 else {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);} 2328 if (patch->patchconstructop == PCPatchConstruct_Star) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);} 2329 else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);} 2330 else if (patch->patchconstructop == PCPatchConstruct_User) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);} 2331 else {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);} 2332 ierr = PetscViewerASCIIPrintf(viewer, "Solver on patches (all same):\n");CHKERRQ(ierr); 2333 if (patch->solver) { 2334 ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 2335 if (!rank) { 2336 ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr); 2337 ierr = PetscObjectView(patch->solver[0], sviewer);CHKERRQ(ierr); 2338 ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr); 2339 } 2340 ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr); 2341 } else { 2342 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 2343 ierr = PetscViewerASCIIPrintf(viewer, "Solver not yet set.\n");CHKERRQ(ierr); 2344 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2345 } 2346 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2347 PetscFunctionReturn(0); 2348 } 2349 2350 /*MC 2351 PCPATCH - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping 2352 small block additive preconditioners. Block definition is based on topology from 2353 a DM and equation numbering from a PetscSection. 2354 2355 Options Database Keys: 2356 + -pc_patch_cells_view - Views the process local cell numbers for each patch 2357 . -pc_patch_points_view - Views the process local mesh point numbers for each patch 2358 . -pc_patch_g2l_view - Views the map between global dofs and patch local dofs for each patch 2359 . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary 2360 - -pc_patch_sub_mat_view - Views the matrix associated with each patch 2361 2362 Level: intermediate 2363 2364 .seealso: PCType, PCCreate(), PCSetType() 2365 M*/ 2366 PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc) 2367 { 2368 PC_PATCH *patch; 2369 PetscErrorCode ierr; 2370 2371 PetscFunctionBegin; 2372 ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr); 2373 2374 if (patch->subspaces_to_exclude) { 2375 PetscHSetIDestroy(&patch->subspaces_to_exclude); 2376 } 2377 PetscHSetICreate(&patch->subspaces_to_exclude); 2378 2379 patch->classname = "pc"; 2380 2381 /* Set some defaults */ 2382 patch->combined = PETSC_FALSE; 2383 patch->save_operators = PETSC_TRUE; 2384 patch->local_composition_type = PC_COMPOSITE_ADDITIVE; 2385 patch->partition_of_unity = PETSC_FALSE; 2386 patch->codim = -1; 2387 patch->dim = -1; 2388 patch->vankadim = -1; 2389 patch->ignoredim = -1; 2390 patch->patchconstructop = PCPatchConstruct_Star; 2391 patch->symmetrise_sweep = PETSC_FALSE; 2392 patch->npatch = 0; 2393 patch->userIS = NULL; 2394 patch->optionsSet = PETSC_FALSE; 2395 patch->iterationSet = NULL; 2396 patch->user_patches = PETSC_FALSE; 2397 ierr = PetscStrallocpy(MATDENSE, (char **) &patch->sub_mat_type);CHKERRQ(ierr); 2398 patch->viewPatches = PETSC_FALSE; 2399 patch->viewCells = PETSC_FALSE; 2400 patch->viewPoints = PETSC_FALSE; 2401 patch->viewSection = PETSC_FALSE; 2402 patch->viewMatrix = PETSC_FALSE; 2403 patch->setupsolver = PCSetUp_PATCH_Linear; 2404 patch->applysolver = PCApply_PATCH_Linear; 2405 patch->resetsolver = PCReset_PATCH_Linear; 2406 patch->destroysolver = PCDestroy_PATCH_Linear; 2407 patch->updatemultiplicative = PCUpdateMultiplicative_PATCH_Linear; 2408 2409 pc->data = (void *) patch; 2410 pc->ops->apply = PCApply_PATCH; 2411 pc->ops->applytranspose = 0; /* PCApplyTranspose_PATCH; */ 2412 pc->ops->setup = PCSetUp_PATCH; 2413 pc->ops->reset = PCReset_PATCH; 2414 pc->ops->destroy = PCDestroy_PATCH; 2415 pc->ops->setfromoptions = PCSetFromOptions_PATCH; 2416 pc->ops->setuponblocks = PCSetUpOnBlocks_PATCH; 2417 pc->ops->view = PCView_PATCH; 2418 pc->ops->applyrichardson = 0; 2419 2420 PetscFunctionReturn(0); 2421 } 2422