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