xref: /petsc/src/ksp/pc/impls/patch/pcpatch.c (revision 92d509846dbdda9cf52b8978c5a86e667c1b6395)
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 PCPatchSetSubMatType(PC pc, MatType sub_mat_type)
287 {
288   PC_PATCH      *patch = (PC_PATCH *) pc->data;
289   PetscErrorCode ierr;
290 
291   PetscFunctionBegin;
292   if (patch->sub_mat_type) {ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);}
293   ierr = PetscStrallocpy(sub_mat_type, (char **) &patch->sub_mat_type);CHKERRQ(ierr);
294   PetscFunctionReturn(0);
295 }
296 
297 /* TODO: Docs */
298 PetscErrorCode PCPatchGetSubMatType(PC pc, MatType *sub_mat_type)
299 {
300   PC_PATCH *patch = (PC_PATCH *) pc->data;
301   PetscFunctionBegin;
302   *sub_mat_type = patch->sub_mat_type;
303   PetscFunctionReturn(0);
304 }
305 
306 /* TODO: Docs */
307 PetscErrorCode PCPatchSetCellNumbering(PC pc, PetscSection cellNumbering)
308 {
309   PC_PATCH      *patch = (PC_PATCH *) pc->data;
310   PetscErrorCode ierr;
311 
312   PetscFunctionBegin;
313   patch->cellNumbering = cellNumbering;
314   ierr = PetscObjectReference((PetscObject) cellNumbering);CHKERRQ(ierr);
315   PetscFunctionReturn(0);
316 }
317 
318 /* TODO: Docs */
319 PetscErrorCode PCPatchGetCellNumbering(PC pc, PetscSection *cellNumbering)
320 {
321   PC_PATCH *patch = (PC_PATCH *) pc->data;
322   PetscFunctionBegin;
323   *cellNumbering = patch->cellNumbering;
324   PetscFunctionReturn(0);
325 }
326 
327 /* TODO: Docs */
328 PetscErrorCode PCPatchSetConstructType(PC pc, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx)
329 {
330   PC_PATCH *patch = (PC_PATCH *) pc->data;
331 
332   PetscFunctionBegin;
333   patch->ctype = ctype;
334   switch (ctype) {
335   case PC_PATCH_STAR:
336     patch->user_patches     = PETSC_FALSE;
337     patch->patchconstructop = PCPatchConstruct_Star;
338     break;
339   case PC_PATCH_VANKA:
340     patch->user_patches     = PETSC_FALSE;
341     patch->patchconstructop = PCPatchConstruct_Vanka;
342     break;
343   case PC_PATCH_USER:
344   case PC_PATCH_PYTHON:
345     patch->user_patches     = PETSC_TRUE;
346     patch->patchconstructop = PCPatchConstruct_User;
347     if (func) {
348       patch->userpatchconstructionop = func;
349       patch->userpatchconstructctx   = ctx;
350     }
351     break;
352   default:
353     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype);
354   }
355   PetscFunctionReturn(0);
356 }
357 
358 /* TODO: Docs */
359 PetscErrorCode PCPatchGetConstructType(PC pc, PCPatchConstructType *ctype, PetscErrorCode (**func)(PC, PetscInt *, IS **, IS *, void *), void **ctx)
360 {
361   PC_PATCH *patch = (PC_PATCH *) pc->data;
362 
363   PetscFunctionBegin;
364   *ctype = patch->ctype;
365   switch (patch->ctype) {
366   case PC_PATCH_STAR:
367   case PC_PATCH_VANKA:
368     break;
369   case PC_PATCH_USER:
370   case PC_PATCH_PYTHON:
371     *func = patch->userpatchconstructionop;
372     *ctx  = patch->userpatchconstructctx;
373     break;
374   default:
375     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype);
376   }
377   PetscFunctionReturn(0);
378 }
379 
380 /* TODO: Docs */
381 PetscErrorCode PCPatchSetDiscretisationInfo(PC pc, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap,
382                                             const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes)
383 {
384   PC_PATCH      *patch = (PC_PATCH *) pc->data;
385   DM             dm;
386   PetscSF       *sfs;
387   PetscInt       cStart, cEnd, i, j;
388   PetscErrorCode ierr;
389 
390   PetscFunctionBegin;
391   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
392   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
393   ierr = PetscMalloc1(nsubspaces, &sfs);CHKERRQ(ierr);
394   ierr = PetscMalloc1(nsubspaces, &patch->dofSection);CHKERRQ(ierr);
395   ierr = PetscMalloc1(nsubspaces, &patch->bs);CHKERRQ(ierr);
396   ierr = PetscMalloc1(nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr);
397   ierr = PetscMalloc1(nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr);
398   ierr = PetscMalloc1(nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr);
399 
400   patch->nsubspaces       = nsubspaces;
401   patch->totalDofsPerCell = 0;
402   for (i = 0; i < nsubspaces; ++i) {
403     ierr = DMGetDefaultSection(dms[i], &patch->dofSection[i]);CHKERRQ(ierr);
404     ierr = PetscObjectReference((PetscObject) patch->dofSection[i]);CHKERRQ(ierr);
405     ierr = DMGetDefaultSF(dms[i], &sfs[i]);CHKERRQ(ierr);
406     patch->bs[i]              = bs[i];
407     patch->nodesPerCell[i]    = nodesPerCell[i];
408     patch->totalDofsPerCell  += nodesPerCell[i]*bs[i];
409     ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i]*bs[i], &patch->cellNodeMap[i]);CHKERRQ(ierr);
410     for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]*bs[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j];
411     patch->subspaceOffsets[i] = subspaceOffsets[i];
412   }
413   ierr = PCPatchCreateDefaultSF_Private(pc, nsubspaces, sfs, patch->bs);CHKERRQ(ierr);
414   ierr = PetscFree(sfs);CHKERRQ(ierr);
415 
416   patch->subspaceOffsets[nsubspaces] = subspaceOffsets[nsubspaces];
417   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr);
418   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr);
419   PetscFunctionReturn(0);
420 }
421 
422 /* TODO: Docs */
423 PetscErrorCode PCPatchSetDiscretisationInfoCombined(PC pc, DM dm, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes)
424 {
425   PC_PATCH      *patch = (PC_PATCH *) pc->data;
426   PetscInt       cStart, cEnd, i, j;
427   PetscErrorCode ierr;
428 
429   PetscFunctionBegin;
430   patch->combined = PETSC_TRUE;
431   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
432   ierr = DMGetNumFields(dm, &patch->nsubspaces);CHKERRQ(ierr);
433   ierr = PetscCalloc1(patch->nsubspaces, &patch->dofSection);CHKERRQ(ierr);
434   ierr = PetscMalloc1(patch->nsubspaces, &patch->bs);CHKERRQ(ierr);
435   ierr = PetscMalloc1(patch->nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr);
436   ierr = PetscMalloc1(patch->nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr);
437   ierr = PetscCalloc1(patch->nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr);
438   ierr = DMGetDefaultSection(dm, &patch->dofSection[0]);CHKERRQ(ierr);
439   ierr = PetscObjectReference((PetscObject) patch->dofSection[0]);CHKERRQ(ierr);
440   ierr = PetscSectionGetStorageSize(patch->dofSection[0], &patch->subspaceOffsets[patch->nsubspaces]);CHKERRQ(ierr);
441   patch->totalDofsPerCell = 0;
442   for (i = 0; i < patch->nsubspaces; ++i) {
443     patch->bs[i]             = 1;
444     patch->nodesPerCell[i]   = nodesPerCell[i];
445     patch->totalDofsPerCell += nodesPerCell[i];
446     ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr);
447     for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j];
448   }
449   ierr = DMGetDefaultSF(dm, &patch->defaultSF);CHKERRQ(ierr);
450   ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr);
451   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr);
452   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr);
453   PetscFunctionReturn(0);
454 }
455 
456 /*@C
457 
458   PCPatchSetComputeFunction - Set the callback used to compute patch residuals
459 
460   Input Parameters:
461 + pc   - The PC
462 . func - The callback
463 - ctx  - The user context
464 
465   Level: advanced
466 
467   Note:
468   The callback has signature:
469 +  usercomputef(pc, point, x, f, cellIS, n, u, ctx)
470 +  pc     - The PC
471 +  point  - The point
472 +  x      - The input solution (not used in linear problems)
473 +  f      - The patch residual vector
474 +  cellIS - An array of the cell numbers
475 +  n      - The size of g2l
476 +  g2l    - The global to local dof translation table
477 +  ctx    - The user context
478   and can assume that the matrix entries have been set to zero before the call.
479 
480 .seealso: PCPatchSetComputeOperator(), PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo()
481 @*/
482 PetscErrorCode PCPatchSetComputeFunction(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, void *), void *ctx)
483 {
484   PC_PATCH *patch = (PC_PATCH *) pc->data;
485 
486   PetscFunctionBegin;
487   patch->usercomputef    = func;
488   patch->usercomputefctx = ctx;
489   PetscFunctionReturn(0);
490 }
491 
492 /*@C
493 
494   PCPatchSetComputeOperator - Set the callback used to compute patch matrices
495 
496   Input Parameters:
497 + pc   - The PC
498 . func - The callback
499 - ctx  - The user context
500 
501   Level: advanced
502 
503   Note:
504   The callback has signature:
505 +  usercomputeop(pc, point, x, mat, cellIS, n, u, ctx)
506 +  pc     - The PC
507 +  point  - The point
508 +  x      - The input solution (not used in linear problems)
509 +  mat    - The patch matrix
510 +  cellIS - An array of the cell numbers
511 +  n      - The size of g2l
512 +  g2l    - The global to local dof translation table
513 +  ctx    - The user context
514   and can assume that the matrix entries have been set to zero before the call.
515 
516 .seealso: PCPatchGetComputeOperator(), PCPatchSetComputeFunction(), PCPatchSetDiscretisationInfo()
517 @*/
518 PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, void *), void *ctx)
519 {
520   PC_PATCH *patch = (PC_PATCH *) pc->data;
521 
522   PetscFunctionBegin;
523   patch->usercomputeop    = func;
524   patch->usercomputeopctx = ctx;
525   PetscFunctionReturn(0);
526 }
527 
528 /* On entry, ht contains the topological entities whose dofs we are responsible for solving for;
529    on exit, cht contains all the topological entities we need to compute their residuals.
530    In full generality this should incorporate knowledge of the sparsity pattern of the matrix;
531    here we assume a standard FE sparsity pattern.*/
532 /* TODO: Use DMPlexGetAdjacency() */
533 static PetscErrorCode PCPatchCompleteCellPatch(PC pc, PetscHSetI ht, PetscHSetI cht)
534 {
535   DM             dm;
536   PetscHashIter  hi;
537   PetscInt       point;
538   PetscInt      *star = NULL, *closure = NULL;
539   PetscInt       ignoredim, iStart = 0, iEnd = -1, starSize, closureSize, si, ci;
540   PetscErrorCode ierr;
541 
542   PetscFunctionBegin;
543   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
544   ierr = PCPatchGetIgnoreDim(pc, &ignoredim);CHKERRQ(ierr);
545   if (ignoredim >= 0) {ierr = DMPlexGetDepthStratum(dm, ignoredim, &iStart, &iEnd);CHKERRQ(ierr);}
546   ierr = PetscHSetIClear(cht);CHKERRQ(ierr);
547   PetscHashIterBegin(ht, hi);
548   while (!PetscHashIterAtEnd(ht, hi)) {
549 
550     PetscHashIterGetKey(ht, hi, point);
551     PetscHashIterNext(ht, hi);
552 
553     /* Loop over all the cells that this point connects to */
554     ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
555     for (si = 0; si < starSize*2; si += 2) {
556       const PetscInt ownedpoint = star[si];
557       /* TODO Check for point in cht before running through closure again */
558       /* now loop over all entities in the closure of that cell */
559       ierr = DMPlexGetTransitiveClosure(dm, ownedpoint, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
560       for (ci = 0; ci < closureSize*2; ci += 2) {
561         const PetscInt seenpoint = closure[ci];
562         if (ignoredim >= 0 && seenpoint >= iStart && seenpoint < iEnd) continue;
563         ierr = PetscHSetIAdd(cht, seenpoint);CHKERRQ(ierr);
564       }
565     }
566   }
567   ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr);
568   ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_FALSE, NULL, &star);CHKERRQ(ierr);
569   PetscFunctionReturn(0);
570 }
571 
572 static PetscErrorCode PCPatchGetGlobalDofs(PC pc, PetscSection dofSection[], PetscInt f, PetscBool combined, PetscInt p, PetscInt *dof, PetscInt *off)
573 {
574   PetscErrorCode ierr;
575 
576   PetscFunctionBegin;
577   if (combined) {
578     if (f < 0) {
579       if (dof) {ierr = PetscSectionGetDof(dofSection[0], p, dof);CHKERRQ(ierr);}
580       if (off) {ierr = PetscSectionGetOffset(dofSection[0], p, off);CHKERRQ(ierr);}
581     } else {
582       if (dof) {ierr = PetscSectionGetFieldDof(dofSection[0], p, f, dof);CHKERRQ(ierr);}
583       if (off) {ierr = PetscSectionGetFieldOffset(dofSection[0], p, f, off);CHKERRQ(ierr);}
584     }
585   } else {
586     if (f < 0) {
587       PC_PATCH *patch = (PC_PATCH *) pc->data;
588       PetscInt  fdof, g;
589 
590       if (dof) {
591         *dof = 0;
592         for (g = 0; g < patch->nsubspaces; ++g) {
593           ierr = PetscSectionGetDof(dofSection[g], p, &fdof);CHKERRQ(ierr);
594           *dof += fdof;
595         }
596       }
597       if (off) {
598         *off = 0;
599         for (g = 0; g < patch->nsubspaces; ++g) {
600           ierr = PetscSectionGetOffset(dofSection[g], p, &fdof);CHKERRQ(ierr);
601           *off += fdof;
602         }
603       }
604     } else {
605       if (dof) {ierr = PetscSectionGetDof(dofSection[f], p, dof);CHKERRQ(ierr);}
606       if (off) {ierr = PetscSectionGetOffset(dofSection[f], p, off);CHKERRQ(ierr);}
607     }
608   }
609   PetscFunctionReturn(0);
610 }
611 
612 /* Given a hash table with a set of topological entities (pts), compute the degrees of
613    freedom in global concatenated numbering on those entities.
614    For Vanka smoothing, this needs to do something special: ignore dofs of the
615    constraint subspace on entities that aren't the base entity we're building the patch
616    around. */
617 static PetscErrorCode PCPatchGetPointDofs(PC pc, PetscHSetI pts, PetscHSetI dofs, PetscInt base, PetscHSetI* subspaces_to_exclude)
618 {
619   PC_PATCH      *patch = (PC_PATCH *) pc->data;
620   PetscHashIter  hi;
621   PetscInt       ldof, loff;
622   PetscInt       k, p;
623   PetscErrorCode ierr;
624 
625   PetscFunctionBegin;
626   ierr = PetscHSetIClear(dofs);CHKERRQ(ierr);
627   for (k = 0; k < patch->nsubspaces; ++k) {
628     PetscInt subspaceOffset = patch->subspaceOffsets[k];
629     PetscInt bs             = patch->bs[k];
630     PetscInt j, l;
631 
632     if (subspaces_to_exclude != NULL) {
633       PetscBool should_exclude_k = PETSC_FALSE;
634       PetscHSetIHas(*subspaces_to_exclude, k, &should_exclude_k);
635       if (should_exclude_k) {
636         /* only get this subspace dofs at the base entity, not any others */
637         ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, base, &ldof, &loff);CHKERRQ(ierr);
638         if (0 == ldof) continue;
639         for (j = loff; j < ldof + loff; ++j) {
640           for (l = 0; l < bs; ++l) {
641             PetscInt dof = bs*j + l + subspaceOffset;
642             ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr);
643           }
644         }
645         continue; /* skip the other dofs of this subspace */
646       }
647     }
648 
649     PetscHashIterBegin(pts, hi);
650     while (!PetscHashIterAtEnd(pts, hi)) {
651       PetscHashIterGetKey(pts, hi, p);
652       PetscHashIterNext(pts, hi);
653       ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, p, &ldof, &loff);CHKERRQ(ierr);
654       if (0 == ldof) continue;
655       for (j = loff; j < ldof + loff; ++j) {
656         for (l = 0; l < bs; ++l) {
657           PetscInt dof = bs*j + l + subspaceOffset;
658           ierr = PetscHSetIAdd(dofs, dof);CHKERRQ(ierr);
659         }
660       }
661     }
662   }
663   PetscFunctionReturn(0);
664 }
665 
666 /* Given two hash tables A and B, compute the keys in B that are not in A, and put them in C */
667 static PetscErrorCode PCPatchComputeSetDifference_Private(PetscHSetI A, PetscHSetI B, PetscHSetI C)
668 {
669   PetscHashIter  hi;
670   PetscInt       key;
671   PetscBool      flg;
672   PetscErrorCode ierr;
673 
674   PetscFunctionBegin;
675   ierr = PetscHSetIClear(C);CHKERRQ(ierr);
676   PetscHashIterBegin(B, hi);
677   while (!PetscHashIterAtEnd(B, hi)) {
678     PetscHashIterGetKey(B, hi, key);
679     PetscHashIterNext(B, hi);
680     ierr = PetscHSetIHas(A, key, &flg);CHKERRQ(ierr);
681     if (!flg) {ierr = PetscHSetIAdd(C, key);CHKERRQ(ierr);}
682   }
683   PetscFunctionReturn(0);
684 }
685 
686 /*
687  * PCPatchCreateCellPatches - create patches.
688  *
689  * Input Parameters:
690  * + dm - The DMPlex object defining the mesh
691  *
692  * Output Parameters:
693  * + cellCounts  - Section with counts of cells around each vertex
694  * . cells       - IS of the cell point indices of cells in each patch
695  * . pointCounts - Section with counts of cells around each vertex
696  * - point       - IS of the cell point indices of cells in each patch
697  */
698 static PetscErrorCode PCPatchCreateCellPatches(PC pc)
699 {
700   PC_PATCH       *patch = (PC_PATCH *) pc->data;
701   DMLabel         ghost = NULL;
702   DM              dm, plex;
703   PetscHSetI      ht, cht;
704   PetscSection    cellCounts,  pointCounts;
705   PetscInt       *cellsArray, *pointsArray;
706   PetscInt        numCells,    numPoints;
707   const PetscInt *leaves;
708   PetscInt        nleaves, pStart, pEnd, cStart, cEnd, vStart, vEnd, v;
709   PetscBool       isFiredrake;
710   PetscErrorCode  ierr;
711 
712   PetscFunctionBegin;
713   /* Used to keep track of the cells in the patch. */
714   ierr = PetscHSetICreate(&ht);CHKERRQ(ierr);
715   ierr = PetscHSetICreate(&cht);CHKERRQ(ierr);
716 
717   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
718   if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch PC\n");
719   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
720   ierr = DMPlexGetChart(plex, &pStart, &pEnd);CHKERRQ(ierr);
721   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
722 
723   if (patch->user_patches) {
724     ierr = patch->userpatchconstructionop(pc, &patch->npatch, &patch->userIS, &patch->iterationSet, patch->userpatchconstructctx);CHKERRQ(ierr);
725     vStart = 0; vEnd = patch->npatch;
726   } else if (patch->codim < 0) {
727     if (patch->dim < 0) {ierr = DMPlexGetDepthStratum(plex,  0,            &vStart, &vEnd);CHKERRQ(ierr);}
728     else                {ierr = DMPlexGetDepthStratum(plex,  patch->dim,   &vStart, &vEnd);CHKERRQ(ierr);}
729   } else                {ierr = DMPlexGetHeightStratum(plex, patch->codim, &vStart, &vEnd);CHKERRQ(ierr);}
730   patch->npatch = vEnd - vStart;
731 
732   /* These labels mark the owned points.  We only create patches around points that this process owns. */
733   ierr = DMHasLabel(dm, "pyop2_ghost", &isFiredrake);CHKERRQ(ierr);
734   if (isFiredrake) {
735     ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr);
736     ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr);
737   } else {
738     PetscSF sf;
739 
740     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
741     ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr);
742     nleaves = PetscMax(nleaves, 0);
743   }
744 
745   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->cellCounts);CHKERRQ(ierr);
746   ierr = PetscObjectSetName((PetscObject) patch->cellCounts, "Patch Cell Layout");CHKERRQ(ierr);
747   cellCounts = patch->cellCounts;
748   ierr = PetscSectionSetChart(cellCounts, vStart, vEnd);CHKERRQ(ierr);
749   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->pointCounts);CHKERRQ(ierr);
750   ierr = PetscObjectSetName((PetscObject) patch->pointCounts, "Patch Point Layout");CHKERRQ(ierr);
751   pointCounts = patch->pointCounts;
752   ierr = PetscSectionSetChart(pointCounts, vStart, vEnd);CHKERRQ(ierr);
753   /* Count cells and points in the patch surrounding each entity */
754   for (v = vStart; v < vEnd; ++v) {
755     PetscHashIter hi;
756     PetscInt       chtSize, loc = -1;
757     PetscBool      flg;
758 
759     if (!patch->user_patches) {
760       if (ghost) {ierr = DMLabelHasPoint(ghost, v, &flg);CHKERRQ(ierr);}
761       else       {ierr = PetscFindInt(v, nleaves, leaves, &loc);CHKERRQ(ierr); flg = loc >=0 ? PETSC_TRUE : PETSC_FALSE;}
762       /* Not an owned entity, don't make a cell patch. */
763       if (flg) continue;
764     }
765 
766     ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr);
767     ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr);
768     ierr = PetscHSetIGetSize(cht, &chtSize);CHKERRQ(ierr);
769     /* empty patch, continue */
770     if (chtSize == 0) continue;
771 
772     /* safe because size(cht) > 0 from above */
773     PetscHashIterBegin(cht, hi);
774     while (!PetscHashIterAtEnd(cht, hi)) {
775       PetscInt point, pdof;
776 
777       PetscHashIterGetKey(cht, hi, point);
778       ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr);
779       if (pdof)                            {ierr = PetscSectionAddDof(pointCounts, v, 1);CHKERRQ(ierr);}
780       if (point >= cStart && point < cEnd) {ierr = PetscSectionAddDof(cellCounts, v, 1);CHKERRQ(ierr);}
781       PetscHashIterNext(cht, hi);
782     }
783   }
784   if (isFiredrake) {ierr = DMLabelDestroyIndex(ghost);CHKERRQ(ierr);}
785 
786   ierr = PetscSectionSetUp(cellCounts);CHKERRQ(ierr);
787   ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr);
788   ierr = PetscMalloc1(numCells, &cellsArray);CHKERRQ(ierr);
789   ierr = PetscSectionSetUp(pointCounts);CHKERRQ(ierr);
790   ierr = PetscSectionGetStorageSize(pointCounts, &numPoints);CHKERRQ(ierr);
791   ierr = PetscMalloc1(numPoints, &pointsArray);CHKERRQ(ierr);
792 
793   /* Now that we know how much space we need, run through again and actually remember the cells. */
794   for (v = vStart; v < vEnd; v++ ) {
795     PetscHashIter hi;
796     PetscInt       dof, off, cdof, coff, pdof, n = 0, cn = 0;
797 
798     ierr = PetscSectionGetDof(pointCounts, v, &dof);CHKERRQ(ierr);
799     ierr = PetscSectionGetOffset(pointCounts, v, &off);CHKERRQ(ierr);
800     ierr = PetscSectionGetDof(cellCounts, v, &cdof);CHKERRQ(ierr);
801     ierr = PetscSectionGetOffset(cellCounts, v, &coff);CHKERRQ(ierr);
802     if (dof <= 0) continue;
803     ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr);
804     ierr = PCPatchCompleteCellPatch(pc, ht, cht);CHKERRQ(ierr);
805     PetscHashIterBegin(cht, hi);
806     while (!PetscHashIterAtEnd(cht, hi)) {
807       PetscInt point;
808 
809       PetscHashIterGetKey(cht, hi, point);
810       ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, -1, patch->combined, point, &pdof, NULL);CHKERRQ(ierr);
811       if (pdof)                            {pointsArray[off + n++] = point;}
812       if (point >= cStart && point < cEnd) {cellsArray[coff + cn++] = point;}
813       PetscHashIterNext(cht, hi);
814     }
815     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);
816     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);
817   }
818   ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr);
819   ierr = PetscHSetIDestroy(&cht);CHKERRQ(ierr);
820   ierr = DMDestroy(&plex);CHKERRQ(ierr);
821 
822   ierr = ISCreateGeneral(PETSC_COMM_SELF, numCells,  cellsArray,  PETSC_OWN_POINTER, &patch->cells);CHKERRQ(ierr);
823   ierr = PetscObjectSetName((PetscObject) patch->cells,  "Patch Cells");CHKERRQ(ierr);
824   if (patch->viewCells) {
825     ierr = ObjectView((PetscObject) patch->cellCounts, patch->viewerCells, patch->formatCells);CHKERRQ(ierr);
826     ierr = ObjectView((PetscObject) patch->cells,      patch->viewerCells, patch->formatCells);CHKERRQ(ierr);
827   }
828   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints, pointsArray, PETSC_OWN_POINTER, &patch->points);CHKERRQ(ierr);
829   ierr = PetscObjectSetName((PetscObject) patch->points, "Patch Points");CHKERRQ(ierr);
830   if (patch->viewPoints) {
831     ierr = ObjectView((PetscObject) patch->pointCounts, patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr);
832     ierr = ObjectView((PetscObject) patch->points,      patch->viewerPoints, patch->formatPoints);CHKERRQ(ierr);
833   }
834   PetscFunctionReturn(0);
835 }
836 
837 /*
838  * PCPatchCreateCellPatchDiscretisationInfo - Build the dof maps for cell patches
839  *
840  * Input Parameters:
841  * + dm - The DMPlex object defining the mesh
842  * . cellCounts - Section with counts of cells around each vertex
843  * . cells - IS of the cell point indices of cells in each patch
844  * . cellNumbering - Section mapping plex cell points to Firedrake cell indices.
845  * . nodesPerCell - number of nodes per cell.
846  * - cellNodeMap - map from cells to node indices (nodesPerCell * numCells)
847  *
848  * Output Parameters:
849  * + dofs - IS of local dof numbers of each cell in the patch, where local is a patch local numbering
850  * . gtolCounts - Section with counts of dofs per cell patch
851  * - gtol - IS mapping from global dofs to local dofs for each patch.
852  */
853 static PetscErrorCode PCPatchCreateCellPatchDiscretisationInfo(PC pc)
854 {
855   PC_PATCH       *patch           = (PC_PATCH *) pc->data;
856   PetscSection    cellCounts      = patch->cellCounts;
857   PetscSection    pointCounts     = patch->pointCounts;
858   PetscSection    gtolCounts;
859   IS              cells           = patch->cells;
860   IS              points          = patch->points;
861   PetscSection    cellNumbering   = patch->cellNumbering;
862   PetscInt        Nf              = patch->nsubspaces;
863   PetscInt        numCells, numPoints;
864   PetscInt        numDofs;
865   PetscInt        numGlobalDofs;
866   PetscInt        totalDofsPerCell = patch->totalDofsPerCell;
867   PetscInt        vStart, vEnd, v;
868   const PetscInt *cellsArray, *pointsArray;
869   PetscInt       *newCellsArray   = NULL;
870   PetscInt       *dofsArray       = NULL;
871   PetscInt       *offsArray       = NULL;
872   PetscInt       *asmArray        = NULL;
873   PetscInt       *globalDofsArray = NULL;
874   PetscInt        globalIndex     = 0;
875   PetscInt        key             = 0;
876   PetscInt        asmKey          = 0;
877   DM              dm              = NULL;
878   const PetscInt *bcNodes         = NULL;
879   PetscHMapI      ht;
880   PetscHSetI      globalBcs;
881   PetscInt        numBcs;
882   PetscHSetI      ownedpts, seenpts, owneddofs, seendofs, artificialbcs;
883   PetscInt        pStart, pEnd, p, i;
884   PetscErrorCode  ierr;
885 
886   PetscFunctionBegin;
887 
888   ierr = PCGetDM(pc, &dm); CHKERRQ(ierr);
889   /* dofcounts section is cellcounts section * dofPerCell */
890   ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr);
891   ierr = PetscSectionGetStorageSize(patch->pointCounts, &numPoints);CHKERRQ(ierr);
892   numDofs = numCells * totalDofsPerCell;
893   ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr);
894   ierr = PetscMalloc1(numPoints*Nf, &offsArray);CHKERRQ(ierr);
895   ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr);
896   ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr);
897   ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr);
898   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr);
899   gtolCounts = patch->gtolCounts;
900   ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr);
901   ierr = PetscObjectSetName((PetscObject) patch->gtolCounts, "Patch Global Index Section");CHKERRQ(ierr);
902 
903   /* Outside the patch loop, get the dofs that are globally-enforced Dirichlet
904    conditions */
905   ierr = PetscHSetICreate(&globalBcs);CHKERRQ(ierr);
906   ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr);
907   ierr = ISGetSize(patch->ghostBcNodes, &numBcs); CHKERRQ(ierr);
908   for (i = 0; i < numBcs; ++i) {
909     ierr = PetscHSetIAdd(globalBcs, bcNodes[i]);CHKERRQ(ierr); /* these are already in concatenated numbering */
910   }
911   ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes); CHKERRQ(ierr);
912   ierr = ISDestroy(&patch->ghostBcNodes); CHKERRQ(ierr); /* memory optimisation */
913 
914   /* Hash tables for artificial BC construction */
915   ierr = PetscHSetICreate(&ownedpts);CHKERRQ(ierr);
916   ierr = PetscHSetICreate(&seenpts);CHKERRQ(ierr);
917   ierr = PetscHSetICreate(&owneddofs);CHKERRQ(ierr);
918   ierr = PetscHSetICreate(&seendofs);CHKERRQ(ierr);
919   ierr = PetscHSetICreate(&artificialbcs);CHKERRQ(ierr);
920 
921   ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr);
922   ierr = ISGetIndices(points, &pointsArray);CHKERRQ(ierr);
923   ierr = PetscHMapICreate(&ht);CHKERRQ(ierr);
924   for (v = vStart; v < vEnd; ++v) {
925     PetscInt localIndex = 0;
926     PetscInt dof, off, i, j, k, l;
927 
928     ierr = PetscHMapIClear(ht);CHKERRQ(ierr);
929     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
930     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
931     if (dof <= 0) continue;
932 
933     /* Calculate the global numbers of the artificial BC dofs here first */
934     ierr = patch->patchconstructop((void*)patch, dm, v, ownedpts); CHKERRQ(ierr);
935     ierr = PCPatchCompleteCellPatch(pc, ownedpts, seenpts); CHKERRQ(ierr);
936     ierr = PCPatchGetPointDofs(pc, ownedpts, owneddofs, v, &patch->subspaces_to_exclude); CHKERRQ(ierr);
937     ierr = PCPatchGetPointDofs(pc, seenpts, seendofs, v, NULL); CHKERRQ(ierr);
938     ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs); CHKERRQ(ierr);
939     if (patch->viewPatches) {
940       PetscHSetI globalbcdofs;
941       PetscHashIter hi;
942       MPI_Comm comm = PetscObjectComm((PetscObject)pc);
943 
944       ierr = PetscHSetICreate(&globalbcdofs);CHKERRQ(ierr);
945       ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v); CHKERRQ(ierr);
946       PetscHashIterBegin(owneddofs, hi);
947       while (!PetscHashIterAtEnd(owneddofs, hi)) {
948         PetscInt globalDof;
949 
950         PetscHashIterGetKey(owneddofs, hi, globalDof);
951         PetscHashIterNext(owneddofs, hi);
952         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
953       }
954       ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr);
955       ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v); CHKERRQ(ierr);
956       PetscHashIterBegin(seendofs, hi);
957       while (!PetscHashIterAtEnd(seendofs, hi)) {
958         PetscInt globalDof;
959         PetscBool flg;
960 
961         PetscHashIterGetKey(seendofs, hi, globalDof);
962         PetscHashIterNext(seendofs, hi);
963         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
964 
965         ierr = PetscHSetIHas(globalBcs, globalDof, &flg);CHKERRQ(ierr);
966         if (flg) {ierr = PetscHSetIAdd(globalbcdofs, globalDof);CHKERRQ(ierr);}
967       }
968       ierr = PetscSynchronizedPrintf(comm, "\n"); CHKERRQ(ierr);
969       ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v); CHKERRQ(ierr);
970       ierr = PetscHSetIGetSize(globalbcdofs, &numBcs);CHKERRQ(ierr);
971       if (numBcs > 0) {
972         PetscHashIterBegin(globalbcdofs, hi);
973         while (!PetscHashIterAtEnd(globalbcdofs, hi)) {
974           PetscInt globalDof;
975           PetscHashIterGetKey(globalbcdofs, hi, globalDof);
976           PetscHashIterNext(globalbcdofs, hi);
977           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr);
978         }
979       }
980       ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr);
981       ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr);
982       ierr = PetscHSetIGetSize(artificialbcs, &numBcs);CHKERRQ(ierr);
983       if (numBcs > 0) {
984         PetscHashIterBegin(artificialbcs, hi);
985         while (!PetscHashIterAtEnd(artificialbcs, hi)) {
986           PetscInt globalDof;
987           PetscHashIterGetKey(artificialbcs, hi, globalDof);
988           PetscHashIterNext(artificialbcs, hi);
989           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof); CHKERRQ(ierr);
990         }
991       }
992       ierr = PetscSynchronizedPrintf(comm, "\n\n"); CHKERRQ(ierr);
993       ierr = PetscHSetIDestroy(&globalbcdofs);CHKERRQ(ierr);
994     }
995    for (k = 0; k < patch->nsubspaces; ++k) {
996       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
997       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
998       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
999       PetscInt        bs             = patch->bs[k];
1000 
1001       for (i = off; i < off + dof; ++i) {
1002         /* Walk over the cells in this patch. */
1003         const PetscInt c    = cellsArray[i];
1004         PetscInt       cell = c;
1005 
1006         /* TODO Change this to an IS */
1007         if (cellNumbering) {
1008           ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr);
1009           if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c);
1010           ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);
1011         }
1012         newCellsArray[i] = cell;
1013         for (j = 0; j < nodesPerCell; ++j) {
1014           /* For each global dof, map it into contiguous local storage. */
1015           const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset;
1016           /* finally, loop over block size */
1017           for (l = 0; l < bs; ++l) {
1018             PetscInt  localDof;
1019             PetscBool isGlobalBcDof, isArtificialBcDof;
1020 
1021             /* first, check if this is either a globally enforced or locally enforced BC dof */
1022             ierr = PetscHSetIHas(globalBcs, globalDof + l, &isGlobalBcDof);CHKERRQ(ierr);
1023             ierr = PetscHSetIHas(artificialbcs, globalDof + l, &isArtificialBcDof);CHKERRQ(ierr);
1024 
1025             /* if it's either, don't ever give it a local dof number */
1026             if (isGlobalBcDof || isArtificialBcDof) {
1027               dofsArray[globalIndex++] = -1; /* don't use this in assembly in this patch */
1028             } else {
1029               ierr = PetscHMapIGet(ht, globalDof + l, &localDof);CHKERRQ(ierr);
1030               if (localDof == -1) {
1031                 localDof = localIndex++;
1032                 ierr = PetscHMapISet(ht, globalDof + l, localDof);CHKERRQ(ierr);
1033               }
1034               if ( globalIndex >= numDofs ) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs);
1035               /* And store. */
1036               dofsArray[globalIndex++] = localDof;
1037             }
1038           }
1039         }
1040       }
1041     }
1042     /* How many local dofs in this patch? */
1043     ierr = PetscHMapIGetSize(ht, &dof);CHKERRQ(ierr);
1044     ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr);
1045   }
1046   if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex);
1047   ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr);
1048   ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr);
1049   ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr);
1050 
1051   /* Now populate the global to local map.  This could be merged into the above loop if we were willing to deal with reallocs. */
1052   for (v = vStart; v < vEnd; ++v) {
1053     PetscHashIter hi;
1054     PetscInt      dof, off, Np, ooff, i, j, k, l;
1055 
1056     ierr = PetscHMapIClear(ht);CHKERRQ(ierr);
1057     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
1058     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
1059     ierr = PetscSectionGetDof(pointCounts, v, &Np);CHKERRQ(ierr);
1060     ierr = PetscSectionGetOffset(pointCounts, v, &ooff);CHKERRQ(ierr);
1061     if (dof <= 0) continue;
1062 
1063     for (k = 0; k < patch->nsubspaces; ++k) {
1064       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
1065       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
1066       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
1067       PetscInt        bs             = patch->bs[k];
1068       PetscInt        goff;
1069 
1070       for (i = off; i < off + dof; ++i) {
1071         /* Reconstruct mapping of global-to-local on this patch. */
1072         const PetscInt c    = cellsArray[i];
1073         PetscInt       cell = c;
1074 
1075         if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);}
1076         for (j = 0; j < nodesPerCell; ++j) {
1077           for (l = 0; l < bs; ++l) {
1078             const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset;
1079             const PetscInt localDof  = dofsArray[key++];
1080 
1081             if (localDof >= 0) {ierr = PetscHMapISet(ht, globalDof, localDof);CHKERRQ(ierr);}
1082           }
1083         }
1084       }
1085 
1086       /* Shove it in the output data structure. */
1087       ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr);
1088       PetscHashIterBegin(ht, hi);
1089       while (!PetscHashIterAtEnd(ht, hi)) {
1090         PetscInt globalDof, localDof;
1091 
1092         PetscHashIterGetKey(ht, hi, globalDof);
1093         PetscHashIterGetVal(ht, hi, localDof);
1094         if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof;
1095         PetscHashIterNext(ht, hi);
1096       }
1097 
1098       for (p = 0; p < Np; ++p) {
1099         const PetscInt point = pointsArray[ooff + p];
1100         PetscInt       globalDof, localDof;
1101 
1102         ierr = PCPatchGetGlobalDofs(pc, patch->dofSection, k, patch->combined, point, NULL, &globalDof);CHKERRQ(ierr);
1103         ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr);
1104         offsArray[(ooff + p)*Nf + k] = localDof;
1105       }
1106     }
1107 
1108     ierr = PetscHSetIDestroy(&globalBcs);CHKERRQ(ierr);
1109     ierr = PetscHSetIDestroy(&ownedpts);CHKERRQ(ierr);
1110     ierr = PetscHSetIDestroy(&seenpts);CHKERRQ(ierr);
1111     ierr = PetscHSetIDestroy(&owneddofs);CHKERRQ(ierr);
1112     ierr = PetscHSetIDestroy(&seendofs);CHKERRQ(ierr);
1113     ierr = PetscHSetIDestroy(&artificialbcs);CHKERRQ(ierr);
1114 
1115       /* At this point, we have a hash table ht built that maps globalDof -> localDof.
1116      We need to create the dof table laid out cellwise first, then by subspace,
1117      as the assembler assembles cell-wise and we need to stuff the different
1118      contributions of the different function spaces to the right places. So we loop
1119      over cells, then over subspaces. */
1120     if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */
1121       for (i = off; i < off + dof; ++i) {
1122         const PetscInt c    = cellsArray[i];
1123         PetscInt       cell = c;
1124 
1125         if (cellNumbering) {ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);}
1126         for (k = 0; k < patch->nsubspaces; ++k) {
1127           const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
1128           PetscInt        nodesPerCell   = patch->nodesPerCell[k];
1129           PetscInt        subspaceOffset = patch->subspaceOffsets[k];
1130           PetscInt        bs             = patch->bs[k];
1131 
1132           for (j = 0; j < nodesPerCell; ++j) {
1133             for (l = 0; l < bs; ++l) {
1134               const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + l + subspaceOffset;
1135               PetscInt       localDof;
1136 
1137               ierr = PetscHMapIGet(ht, globalDof, &localDof);CHKERRQ(ierr);
1138               /* If it's not in the hash table, i.e. is a BC dof,
1139                then the PetscHSetIMap above gives -1, which matches
1140                exactly the convention for PETSc's matrix assembly to
1141                ignore the dof. So we don't need to do anything here */
1142               asmArray[asmKey++] = localDof;
1143             }
1144           }
1145         }
1146       }
1147     }
1148   }
1149   if (1 == patch->nsubspaces) {ierr = PetscMemcpy(asmArray, dofsArray, numDofs * sizeof(PetscInt));CHKERRQ(ierr);}
1150 
1151   ierr = PetscHMapIDestroy(&ht);CHKERRQ(ierr);
1152   ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr);
1153   ierr = ISRestoreIndices(points, &pointsArray);CHKERRQ(ierr);
1154   ierr = PetscFree(dofsArray);CHKERRQ(ierr);
1155   /* Create placeholder section for map from points to patch dofs */
1156   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->patchSection);CHKERRQ(ierr);
1157   ierr = PetscSectionSetNumFields(patch->patchSection, patch->nsubspaces);CHKERRQ(ierr);
1158   if (patch->combined) {
1159     PetscInt numFields;
1160     ierr = PetscSectionGetNumFields(patch->dofSection[0], &numFields);CHKERRQ(ierr);
1161     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);
1162     ierr = PetscSectionGetChart(patch->dofSection[0], &pStart, &pEnd);CHKERRQ(ierr);
1163     ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr);
1164     for (p = pStart; p < pEnd; ++p) {
1165       PetscInt dof, fdof, f;
1166 
1167       ierr = PetscSectionGetDof(patch->dofSection[0], p, &dof);CHKERRQ(ierr);
1168       ierr = PetscSectionSetDof(patch->patchSection, p, dof);CHKERRQ(ierr);
1169       for (f = 0; f < patch->nsubspaces; ++f) {
1170         ierr = PetscSectionGetFieldDof(patch->dofSection[0], p, f, &fdof);CHKERRQ(ierr);
1171         ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr);
1172       }
1173     }
1174   } else {
1175     PetscInt pStartf, pEndf, f;
1176     pStart = PETSC_MAX_INT;
1177     pEnd = PETSC_MIN_INT;
1178     for (f = 0; f < patch->nsubspaces; ++f) {
1179       ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr);
1180       pStart = PetscMin(pStart, pStartf);
1181       pEnd = PetscMax(pEnd, pEndf);
1182     }
1183     ierr = PetscSectionSetChart(patch->patchSection, pStart, pEnd);CHKERRQ(ierr);
1184     for (f = 0; f < patch->nsubspaces; ++f) {
1185       ierr = PetscSectionGetChart(patch->dofSection[f], &pStartf, &pEndf);CHKERRQ(ierr);
1186       for (p = pStartf; p < pEndf; ++p) {
1187         PetscInt fdof;
1188         ierr = PetscSectionGetDof(patch->dofSection[f], p, &fdof);CHKERRQ(ierr);
1189         ierr = PetscSectionAddDof(patch->patchSection, p, fdof);CHKERRQ(ierr);
1190         ierr = PetscSectionSetFieldDof(patch->patchSection, p, f, fdof);CHKERRQ(ierr);
1191       }
1192     }
1193   }
1194   ierr = PetscSectionSetUp(patch->patchSection);CHKERRQ(ierr);
1195   ierr = PetscSectionSetUseFieldOffsets(patch->patchSection, PETSC_TRUE);CHKERRQ(ierr);
1196   /* Replace cell indices with firedrake-numbered ones. */
1197   ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr);
1198   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr);
1199   ierr = PetscObjectSetName((PetscObject) patch->gtol, "Global Indices");CHKERRQ(ierr);
1200   ierr = PetscSectionViewFromOptions(patch->gtolCounts, (PetscObject) pc, "-pc_patch_g2l_view");CHKERRQ(ierr);
1201   ierr = ISViewFromOptions(patch->gtol, (PetscObject) pc, "-pc_patch_g2l_view");CHKERRQ(ierr);
1202   ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr);
1203   ierr = ISCreateGeneral(PETSC_COMM_SELF, numPoints*Nf, offsArray, PETSC_OWN_POINTER, &patch->offs);CHKERRQ(ierr);
1204   PetscFunctionReturn(0);
1205 }
1206 
1207 static PetscErrorCode PCPatchZeroFillMatrix_Private(Mat mat, const PetscInt ncell, const PetscInt ndof, const PetscInt *dof)
1208 {
1209   PetscScalar    *values = NULL;
1210   PetscInt        rows, c, i;
1211   PetscErrorCode  ierr;
1212 
1213   PetscFunctionBegin;
1214   ierr = PetscCalloc1(ndof*ndof, &values);CHKERRQ(ierr);
1215   for (c = 0; c < ncell; ++c) {
1216     const PetscInt *idx = &dof[ndof*c];
1217     ierr = MatSetValues(mat, ndof, idx, ndof, idx, values, INSERT_VALUES);CHKERRQ(ierr);
1218   }
1219   ierr = MatGetLocalSize(mat, &rows, NULL);CHKERRQ(ierr);
1220   for (i = 0; i < rows; ++i) {
1221     ierr = MatSetValues(mat, 1, &i, 1, &i, values, INSERT_VALUES);CHKERRQ(ierr);
1222   }
1223   ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1224   ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1225   ierr = PetscFree(values);CHKERRQ(ierr);
1226   PetscFunctionReturn(0);
1227 }
1228 
1229 static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat)
1230 {
1231   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1232   Vec            x, y;
1233   PetscBool      flg;
1234   PetscInt       csize, rsize;
1235   const char    *prefix = NULL;
1236   PetscErrorCode ierr;
1237 
1238   PetscFunctionBegin;
1239   x = patch->patchX[point];
1240   y = patch->patchY[point];
1241   ierr = VecGetSize(x, &csize);CHKERRQ(ierr);
1242   ierr = VecGetSize(y, &rsize);CHKERRQ(ierr);
1243   ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr);
1244   ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
1245   ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr);
1246   ierr = MatAppendOptionsPrefix(*mat, "pc_patch_sub_");CHKERRQ(ierr);
1247   if (patch->sub_mat_type)       {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);}
1248   else if (!patch->sub_mat_type) {ierr = MatSetType(*mat, MATDENSE);CHKERRQ(ierr);}
1249   ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr);
1250   ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr);
1251   if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);}
1252   /* Sparse patch matrices */
1253   if (!flg) {
1254     PetscBT         bt;
1255     PetscInt       *dnnz      = NULL;
1256     const PetscInt *dofsArray = NULL;
1257     PetscInt        pStart, pEnd, ncell, offset, c, i, j;
1258 
1259     ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1260     ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1261     point += pStart;
1262     if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1263     ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1264     ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1265     ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr);
1266     ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
1267     /* XXX: This uses N^2 bits to store the sparsity pattern on a
1268      * patch.  This is probably OK if the patches are not too big,
1269      * but could use quite a bit of memory for planes in 3D.
1270      * Should we switch based on the value of rsize to a
1271      * hash-table (slower, but more memory efficient) approach? */
1272     ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr);
1273     for (c = 0; c < ncell; ++c) {
1274       const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell;
1275       for (i = 0; i < patch->totalDofsPerCell; ++i) {
1276         const PetscInt row = idx[i];
1277         if (row < 0) continue;
1278         for (j = 0; j < patch->totalDofsPerCell; ++j) {
1279           const PetscInt col = idx[j];
1280           const PetscInt key = row*rsize + col;
1281           if (col < 0) continue;
1282           if (!PetscBTLookupSet(bt, key)) ++dnnz[row];
1283         }
1284       }
1285     }
1286     ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
1287     ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr);
1288     ierr = PetscFree(dnnz);CHKERRQ(ierr);
1289     ierr = PCPatchZeroFillMatrix_Private(*mat, ncell, patch->totalDofsPerCell, &dofsArray[offset*patch->totalDofsPerCell]);CHKERRQ(ierr);
1290     ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
1291     ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1292   }
1293   ierr = MatSetUp(*mat);CHKERRQ(ierr);
1294   PetscFunctionReturn(0);
1295 }
1296 
1297 static PetscErrorCode PCPatchComputeFunction_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Vec F, IS cellIS, PetscInt n, const PetscInt *l2p, void *ctx)
1298 {
1299   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1300   DM              dm;
1301   PetscSection    s;
1302   const PetscInt *parray, *oarray;
1303   PetscInt        Nf = patch->nsubspaces, Np, poff, p, f;
1304   PetscErrorCode  ierr;
1305 
1306   PetscFunctionBegin;
1307   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1308   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1309   /* Set offset into patch */
1310   ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr);
1311   ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr);
1312   ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr);
1313   ierr = ISGetIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1314   for (f = 0; f < Nf; ++f) {
1315     for (p = 0; p < Np; ++p) {
1316       const PetscInt point = parray[poff+p];
1317       PetscInt       dof;
1318 
1319       ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr);
1320       ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);
1321       if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);}
1322       else                        {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);}
1323     }
1324   }
1325   ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr);
1326   ierr = ISRestoreIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1327   if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);}
1328   ierr = DMPlexComputeResidual_Patch_Internal(pc->dm, patch->patchSection, cellIS, 0.0, x, NULL, F, ctx);CHKERRQ(ierr);
1329   PetscFunctionReturn(0);
1330 }
1331 
1332 PetscErrorCode PCPatchComputeFunction_Internal(PC pc, Vec x, Vec F, PetscInt point)
1333 {
1334   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1335   const PetscInt *dofsArray;
1336   const PetscInt *cellsArray;
1337   PetscInt        ncell, offset, pStart, pEnd;
1338   PetscErrorCode  ierr;
1339 
1340   PetscFunctionBegin;
1341   ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1342   if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n");
1343   ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1344   ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1345   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1346 
1347   point += pStart;
1348   if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1349 
1350   ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1351   ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1352   if (ncell <= 0) {
1353     ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1354     PetscFunctionReturn(0);
1355   }
1356   PetscStackPush("PCPatch user callback");
1357   /* Cannot reuse the same IS because the geometry info is being cached in it */
1358   ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr);
1359   ierr = patch->usercomputef(pc, point, x, F, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputefctx);CHKERRQ(ierr);
1360   PetscStackPop;
1361   ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr);
1362   ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1363   ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1364   if (patch->viewMatrix) {
1365     char name[PETSC_MAX_PATH_LEN];
1366 
1367     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch vector for Point %D", point);CHKERRQ(ierr);
1368     ierr = PetscObjectSetName((PetscObject) F, name);CHKERRQ(ierr);
1369     ierr = ObjectView((PetscObject) F, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr);
1370   }
1371   ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1372   PetscFunctionReturn(0);
1373 }
1374 
1375 static PetscErrorCode PCPatchComputeOperator_DMPlex_Private(PC pc, PetscInt patchNum, Vec x, Mat J, IS cellIS, PetscInt n, const PetscInt *l2p, void *ctx)
1376 {
1377   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1378   DM              dm;
1379   PetscSection    s;
1380   const PetscInt *parray, *oarray;
1381   PetscInt        Nf = patch->nsubspaces, Np, poff, p, f;
1382   PetscErrorCode  ierr;
1383 
1384   PetscFunctionBegin;
1385   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1386   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1387   /* Set offset into patch */
1388   ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr);
1389   ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr);
1390   ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr);
1391   ierr = ISGetIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1392   for (f = 0; f < Nf; ++f) {
1393     for (p = 0; p < Np; ++p) {
1394       const PetscInt point = parray[poff+p];
1395       PetscInt       dof;
1396 
1397       ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr);
1398       ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);
1399       if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);}
1400       else                        {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);}
1401     }
1402   }
1403   ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr);
1404   ierr = ISRestoreIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1405   if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);}
1406   /* TODO Shut off MatViewFromOptions() in MatAssemblyEnd() here */
1407   ierr = DMPlexComputeJacobian_Patch_Internal(pc->dm, patch->patchSection, patch->patchSection, cellIS, 0.0, 0.0, x, NULL, J, J, ctx);CHKERRQ(ierr);
1408   PetscFunctionReturn(0);
1409 }
1410 
1411 PetscErrorCode PCPatchComputeOperator_Internal(PC pc, Vec x, Mat mat, PetscInt point)
1412 {
1413   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1414   const PetscInt *dofsArray;
1415   const PetscInt *cellsArray;
1416   PetscInt        ncell, offset, pStart, pEnd;
1417   PetscErrorCode  ierr;
1418 
1419   PetscFunctionBegin;
1420   ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1421   if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n");
1422   ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1423   ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1424   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1425 
1426   point += pStart;
1427   if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1428 
1429   ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1430   ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1431   if (ncell <= 0) {
1432     ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1433     PetscFunctionReturn(0);
1434   }
1435   PetscStackPush("PCPatch user callback");
1436   /* Cannot reuse the same IS because the geometry info is being cached in it */
1437   ierr = ISCreateGeneral(PETSC_COMM_SELF, ncell, cellsArray + offset, PETSC_USE_POINTER, &patch->cellIS);CHKERRQ(ierr);
1438   ierr = patch->usercomputeop(pc, point, x, mat, patch->cellIS, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputeopctx);CHKERRQ(ierr);
1439   PetscStackPop;
1440   ierr = ISDestroy(&patch->cellIS);CHKERRQ(ierr);
1441   ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1442   ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1443   if (patch->viewMatrix) {
1444     char name[PETSC_MAX_PATH_LEN];
1445 
1446     ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN-1, "Patch matrix for Point %D", point);CHKERRQ(ierr);
1447     ierr = PetscObjectSetName((PetscObject) mat, name);CHKERRQ(ierr);
1448     ierr = ObjectView((PetscObject) mat, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr);
1449   }
1450   ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1451   PetscFunctionReturn(0);
1452 }
1453 
1454 static PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat)
1455 {
1456   PC_PATCH          *patch     = (PC_PATCH *) pc->data;
1457   const PetscScalar *xArray    = NULL;
1458   PetscScalar       *yArray    = NULL;
1459   const PetscInt    *gtolArray = NULL;
1460   PetscInt           dof, offset, lidx;
1461   PetscErrorCode     ierr;
1462 
1463   PetscFunctionBeginHot;
1464   ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
1465   ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr);
1466   ierr = VecGetArray(y, &yArray);CHKERRQ(ierr);
1467   ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
1468   ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr);
1469   ierr = ISGetIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1470   if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n");
1471   if (mode == ADD_VALUES    && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n");
1472   for (lidx = 0; lidx < dof; ++lidx) {
1473     const PetscInt gidx = gtolArray[offset+lidx];
1474 
1475     if (mode == INSERT_VALUES) yArray[lidx]  = xArray[gidx]; /* Forward */
1476     else                       yArray[gidx] += xArray[lidx]; /* Reverse */
1477   }
1478   ierr = ISRestoreIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1479   ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr);
1480   ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr);
1481   ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
1482   PetscFunctionReturn(0);
1483 }
1484 
1485 static PetscErrorCode PCSetUp_PATCH_Linear(PC pc)
1486 {
1487   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1488   const char    *prefix;
1489   PetscInt       i;
1490   PetscErrorCode ierr;
1491 
1492   PetscFunctionBegin;
1493   if (!pc->setupcalled) {
1494     ierr = PetscMalloc1(patch->npatch, &patch->solver);CHKERRQ(ierr);
1495     ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
1496     for (i = 0; i < patch->npatch; ++i) {
1497       KSP ksp;
1498       PC  subpc;
1499 
1500       ierr = KSPCreate(PETSC_COMM_SELF, &ksp);CHKERRQ(ierr);
1501       ierr = KSPSetOptionsPrefix(ksp, prefix);CHKERRQ(ierr);
1502       ierr = KSPAppendOptionsPrefix(ksp, "sub_");CHKERRQ(ierr);
1503       ierr = PetscObjectIncrementTabLevel((PetscObject) ksp, (PetscObject) pc, 1);CHKERRQ(ierr);
1504       ierr = KSPGetPC(ksp, &subpc);CHKERRQ(ierr);
1505       ierr = PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1);CHKERRQ(ierr);
1506       ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) ksp);CHKERRQ(ierr);
1507       patch->solver[i] = (PetscObject) ksp;
1508     }
1509   }
1510   if (patch->save_operators) {
1511     for (i = 0; i < patch->npatch; ++i) {
1512       ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr);
1513       ierr = PCPatchComputeOperator_Internal(pc, NULL, patch->mat[i], i);CHKERRQ(ierr);
1514       ierr = KSPSetOperators((KSP) patch->solver[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr);
1515     }
1516   }
1517   PetscFunctionReturn(0);
1518 }
1519 
1520 static PetscErrorCode PCSetUp_PATCH(PC pc)
1521 {
1522   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1523   PetscInt       i;
1524   PetscErrorCode ierr;
1525 
1526   PetscFunctionBegin;
1527   if (!pc->setupcalled) {
1528     PetscInt pStart, pEnd, p;
1529     PetscInt localSize;
1530 
1531     ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
1532 
1533     if (!patch->nsubspaces) {
1534       DM           dm;
1535       PetscDS      prob;
1536       PetscSection s;
1537       PetscInt     cStart, cEnd, c, Nf, f, numGlobalBcs = 0, *globalBcs, *Nb, totNb = 0, **cellDofs;
1538 
1539       ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1540       if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Must set DM for PCPATCH or call PCPatchSetDiscretisationInfo()");
1541       ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1542       ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr);
1543       ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1544       for (p = pStart; p < pEnd; ++p) {
1545         PetscInt cdof;
1546         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1547         numGlobalBcs += cdof;
1548       }
1549       ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1550       ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1551       ierr = PetscMalloc3(Nf, &Nb, Nf, &cellDofs, numGlobalBcs, &globalBcs);CHKERRQ(ierr);
1552       for (f = 0; f < Nf; ++f) {
1553         PetscFE        fe;
1554         PetscDualSpace sp;
1555         PetscInt       cdoff = 0;
1556 
1557         ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
1558         /* ierr = PetscFEGetNumComponents(fe, &Nc[f]);CHKERRQ(ierr); */
1559         ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
1560         ierr = PetscDualSpaceGetDimension(sp, &Nb[f]);CHKERRQ(ierr);
1561         totNb += Nb[f];
1562 
1563         ierr = PetscMalloc1((cEnd-cStart)*Nb[f], &cellDofs[f]);CHKERRQ(ierr);
1564         for (c = cStart; c < cEnd; ++c) {
1565           PetscInt *closure = NULL;
1566           PetscInt  clSize  = 0, cl;
1567 
1568           ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
1569           for (cl = 0; cl < clSize*2; cl += 2) {
1570             const PetscInt p = closure[cl];
1571             PetscInt       fdof, d, foff;
1572 
1573             ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
1574             ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
1575             for (d = 0; d < fdof; ++d, ++cdoff) cellDofs[f][cdoff] = foff + d;
1576           }
1577           ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
1578         }
1579         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]);
1580       }
1581       numGlobalBcs = 0;
1582       for (p = pStart; p < pEnd; ++p) {
1583         const PetscInt *ind;
1584         PetscInt        off, cdof, d;
1585 
1586         ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr);
1587         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1588         ierr = PetscSectionGetConstraintIndices(s, p, &ind);CHKERRQ(ierr);
1589         for (d = 0; d < cdof; ++d) globalBcs[numGlobalBcs++] = off + ind[d];
1590       }
1591 
1592       ierr = PCPatchSetDiscretisationInfoCombined(pc, dm, Nb, (const PetscInt **) cellDofs, numGlobalBcs, globalBcs, numGlobalBcs, globalBcs);CHKERRQ(ierr);
1593       for (f = 0; f < Nf; ++f) {
1594         ierr = PetscFree(cellDofs[f]);CHKERRQ(ierr);
1595       }
1596       ierr = PetscFree3(Nb, cellDofs, globalBcs);CHKERRQ(ierr);
1597       ierr = PCPatchSetComputeFunction(pc, PCPatchComputeFunction_DMPlex_Private, NULL);CHKERRQ(ierr);
1598       ierr = PCPatchSetComputeOperator(pc, PCPatchComputeOperator_DMPlex_Private, NULL);CHKERRQ(ierr);
1599     }
1600 
1601     localSize = patch->subspaceOffsets[patch->nsubspaces];
1602     ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localX);CHKERRQ(ierr);
1603     ierr = VecSetUp(patch->localX);CHKERRQ(ierr);
1604     ierr = VecDuplicate(patch->localX, &patch->localY);CHKERRQ(ierr);
1605     ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr);
1606     ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr);
1607 
1608     /* OK, now build the work vectors */
1609     ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr);
1610     ierr = PetscMalloc1(patch->npatch, &patch->patchX);CHKERRQ(ierr);
1611     ierr = PetscMalloc1(patch->npatch, &patch->patchY);CHKERRQ(ierr);
1612     for (p = pStart; p < pEnd; ++p) {
1613       PetscInt dof;
1614 
1615       ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
1616       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchX[p-pStart]);CHKERRQ(ierr);
1617       ierr = VecSetUp(patch->patchX[p-pStart]);CHKERRQ(ierr);
1618       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchY[p-pStart]);CHKERRQ(ierr);
1619       ierr = VecSetUp(patch->patchY[p-pStart]);CHKERRQ(ierr);
1620     }
1621     if (patch->save_operators) {
1622       ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr);
1623       for (i = 0; i < patch->npatch; ++i) {
1624         ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i]);CHKERRQ(ierr);
1625       }
1626     }
1627     ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
1628 
1629     /* If desired, calculate weights for dof multiplicity */
1630     if (patch->partition_of_unity) {
1631       ierr = VecDuplicate(patch->localX, &patch->dof_weights);CHKERRQ(ierr);
1632       for (i = 0; i < patch->npatch; ++i) {
1633         PetscInt dof;
1634 
1635         ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr);
1636         if (dof <= 0) continue;
1637         ierr = VecSet(patch->patchX[i], 1.0);CHKERRQ(ierr);
1638         ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchX[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
1639       }
1640       ierr = VecReciprocal(patch->dof_weights);CHKERRQ(ierr);
1641     }
1642   }
1643   ierr = (*patch->setupsolver)(pc);CHKERRQ(ierr);
1644   PetscFunctionReturn(0);
1645 }
1646 
1647 static PetscErrorCode PCApply_PATCH_Linear(PC pc, PetscInt i, Vec x, Vec y)
1648 {
1649   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1650   KSP            ksp   = (KSP) patch->solver[i];
1651   PetscErrorCode ierr;
1652 
1653   PetscFunctionBegin;
1654   if (!patch->save_operators) {
1655     Mat mat;
1656 
1657     ierr = PCPatchCreateMatrix_Private(pc, i, &mat);CHKERRQ(ierr);
1658     /* Populate operator here. */
1659     ierr = PCPatchComputeOperator_Internal(pc, NULL, mat, i);CHKERRQ(ierr);
1660     ierr = KSPSetOperators(ksp, mat, mat);CHKERRQ(ierr);
1661     /* Drop reference so the KSPSetOperators below will blow it away. */
1662     ierr = MatDestroy(&mat);CHKERRQ(ierr);
1663   }
1664   ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1665   if (!ksp->setfromoptionscalled) {
1666     ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
1667   }
1668   ierr = KSPSolve(ksp, x, y);CHKERRQ(ierr);
1669   ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1670   if (!patch->save_operators) {
1671     PC pc;
1672     ierr = KSPSetOperators(ksp, NULL, NULL);CHKERRQ(ierr);
1673     ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
1674     /* Destroy PC context too, otherwise the factored matrix hangs around. */
1675     ierr = PCReset(pc);CHKERRQ(ierr);
1676   }
1677   PetscFunctionReturn(0);
1678 }
1679 
1680 static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y)
1681 {
1682   PC_PATCH          *patch    = (PC_PATCH *) pc->data;
1683   const PetscScalar *globalX  = NULL;
1684   PetscScalar       *localX   = NULL;
1685   PetscScalar       *globalY  = NULL;
1686   const PetscInt    *bcNodes  = NULL;
1687   PetscInt           nsweep   = patch->symmetrise_sweep ? 2 : 1;
1688   PetscInt           start[2] = {0, 0};
1689   PetscInt           end[2]   = {-1, -1};
1690   const PetscInt     inc[2]   = {1, -1};
1691   const PetscScalar *localY;
1692   const PetscInt    *iterationSet;
1693   PetscInt           pStart, numBcs, n, sweep, bc, j;
1694   PetscErrorCode     ierr;
1695 
1696   PetscFunctionBegin;
1697   ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
1698   ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr);
1699   /* start, end, inc have 2 entries to manage a second backward sweep if we symmetrize */
1700   end[0]   = patch->npatch;
1701   start[1] = patch->npatch-1;
1702   if (patch->user_patches) {
1703     ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr);
1704     start[1] = end[0] - 1;
1705     ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);
1706   }
1707   /* Scatter from global space into overlapped local spaces */
1708   ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr);
1709   ierr = VecGetArray(patch->localX, &localX);CHKERRQ(ierr);
1710   ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr);
1711   ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr);
1712   ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr);
1713   ierr = VecRestoreArray(patch->localX, &localX);CHKERRQ(ierr);
1714 
1715   ierr = VecSet(patch->localY, 0.0);CHKERRQ(ierr);
1716   ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr);
1717   for (sweep = 0; sweep < nsweep; sweep++) {
1718     for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) {
1719       PetscInt i       = patch->user_patches ? iterationSet[j] : j;
1720       PetscInt start, len;
1721 
1722       ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr);
1723       ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr);
1724       /* TODO: Squash out these guys in the setup as well. */
1725       if (len <= 0) continue;
1726       /* TODO: Do we need different scatters for X and Y? */
1727       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localX, patch->patchX[i], INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
1728       ierr = (*patch->applysolver)(pc, i, patch->patchX[i], patch->patchY[i]);CHKERRQ(ierr);
1729       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchY[i], patch->localY, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
1730     }
1731   }
1732   if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);}
1733   /* XXX: should we do this on the global vector? */
1734   if (patch->partition_of_unity) {
1735     ierr = VecPointwiseMult(patch->localY, patch->localY, patch->dof_weights);CHKERRQ(ierr);
1736   }
1737   /* Now patch->localY contains the solution of the patch solves, so we need to combine them all. */
1738   ierr = VecSet(y, 0.0);CHKERRQ(ierr);
1739   ierr = VecGetArray(y, &globalY);CHKERRQ(ierr);
1740   ierr = VecGetArrayRead(patch->localY, &localY);CHKERRQ(ierr);
1741   ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr);
1742   ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr);
1743   ierr = VecRestoreArrayRead(patch->localY, &localY);CHKERRQ(ierr);
1744 
1745   /* Now we need to send the global BC values through */
1746   ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr);
1747   ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr);
1748   ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
1749   ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr);
1750   for (bc = 0; bc < numBcs; ++bc) {
1751     const PetscInt idx = bcNodes[bc];
1752     if (idx < n) globalY[idx] = globalX[idx];
1753   }
1754 
1755   ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
1756   ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr);
1757   ierr = VecRestoreArray(y, &globalY);CHKERRQ(ierr);
1758 
1759   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
1760   ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
1761   PetscFunctionReturn(0);
1762 }
1763 
1764 static PetscErrorCode PCReset_PATCH_Linear(PC pc)
1765 {
1766   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1767   PetscInt       i;
1768   PetscErrorCode ierr;
1769 
1770   PetscFunctionBegin;
1771   if (patch->solver) {
1772     for (i = 0; i < patch->npatch; ++i) {ierr = KSPReset((KSP) patch->solver[i]);CHKERRQ(ierr);}
1773   }
1774   PetscFunctionReturn(0);
1775 }
1776 
1777 static PetscErrorCode PCReset_PATCH(PC pc)
1778 {
1779   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1780   PetscInt       i;
1781   PetscErrorCode ierr;
1782 
1783   PetscFunctionBegin;
1784   /* TODO: Get rid of all these ifs */
1785   ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr);
1786   ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr);
1787   ierr = PetscSectionDestroy(&patch->pointCounts);CHKERRQ(ierr);
1788   ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr);
1789   ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr);
1790   ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr);
1791   ierr = ISDestroy(&patch->cells);CHKERRQ(ierr);
1792   ierr = ISDestroy(&patch->points);CHKERRQ(ierr);
1793   ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr);
1794   ierr = ISDestroy(&patch->offs);CHKERRQ(ierr);
1795   ierr = PetscSectionDestroy(&patch->patchSection);CHKERRQ(ierr);
1796   ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr);
1797   ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr);
1798 
1799   if (patch->dofSection) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);}
1800   ierr = PetscFree(patch->dofSection);CHKERRQ(ierr);
1801   ierr = PetscFree(patch->bs);CHKERRQ(ierr);
1802   ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr);
1803   if (patch->cellNodeMap) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscFree(patch->cellNodeMap[i]);CHKERRQ(ierr);}
1804   ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr);
1805   ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr);
1806 
1807   ierr = (*patch->resetsolver)(pc);CHKERRQ(ierr);
1808 
1809   if (patch->subspaces_to_exclude) {
1810     PetscHSetIDestroy(&patch->subspaces_to_exclude);
1811   }
1812 
1813   ierr = VecDestroy(&patch->localX);CHKERRQ(ierr);
1814   ierr = VecDestroy(&patch->localY);CHKERRQ(ierr);
1815   if (patch->patchX) {
1816     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchX[i]);CHKERRQ(ierr);}
1817     ierr = PetscFree(patch->patchX);CHKERRQ(ierr);
1818   }
1819   if (patch->patchY) {
1820     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchY[i]);CHKERRQ(ierr);}
1821     ierr = PetscFree(patch->patchY);CHKERRQ(ierr);
1822   }
1823   ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr);
1824   if (patch->patch_dof_weights) {
1825     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);}
1826     ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr);
1827   }
1828   if (patch->mat) {
1829     for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);}
1830     ierr = PetscFree(patch->mat);CHKERRQ(ierr);
1831   }
1832   ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);
1833   if (patch->userIS) {
1834     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);}
1835     ierr = PetscFree(patch->userIS);CHKERRQ(ierr);
1836   }
1837   patch->bs          = 0;
1838   patch->cellNodeMap = NULL;
1839   patch->nsubspaces  = 0;
1840   ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr);
1841 
1842   ierr = PetscViewerDestroy(&patch->viewerSection);CHKERRQ(ierr);
1843   PetscFunctionReturn(0);
1844 }
1845 
1846 static PetscErrorCode PCDestroy_PATCH_Linear(PC pc)
1847 {
1848   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1849   PetscInt       i;
1850   PetscErrorCode ierr;
1851 
1852   PetscFunctionBegin;
1853   if (patch->solver) {
1854     for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy((KSP *) &patch->solver[i]);CHKERRQ(ierr);}
1855     ierr = PetscFree(patch->solver);CHKERRQ(ierr);
1856   }
1857   PetscFunctionReturn(0);
1858 }
1859 
1860 static PetscErrorCode PCDestroy_PATCH(PC pc)
1861 {
1862   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1863   PetscErrorCode ierr;
1864 
1865   PetscFunctionBegin;
1866   ierr = PCReset_PATCH(pc);CHKERRQ(ierr);
1867   ierr = (*patch->destroysolver)(pc);CHKERRQ(ierr);
1868   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1869   PetscFunctionReturn(0);
1870 }
1871 
1872 static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc)
1873 {
1874   PC_PATCH            *patch = (PC_PATCH *) pc->data;
1875   PCPatchConstructType patchConstructionType = PC_PATCH_STAR;
1876   char                 sub_mat_type[PETSC_MAX_PATH_LEN];
1877   const char          *prefix;
1878   PetscBool            flg, dimflg, codimflg;
1879   MPI_Comm             comm;
1880   PetscInt            *ifields, nfields, k;
1881   PetscErrorCode       ierr;
1882 
1883   PetscFunctionBegin;
1884   ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr);
1885   ierr = PetscObjectGetOptionsPrefix((PetscObject) pc, &prefix);CHKERRQ(ierr);
1886   ierr = PetscOptionsHead(PetscOptionsObject, "Vertex-patch Additive Schwarz options");CHKERRQ(ierr);
1887   ierr = PetscOptionsBool("-pc_patch_save_operators",  "Store all patch operators for lifetime of PC?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr);
1888   ierr = PetscOptionsBool("-pc_patch_partition_of_unity", "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr);
1889   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);
1890   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);
1891   if (dimflg && codimflg) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension");CHKERRQ(ierr);
1892   ierr = PetscOptionsEnum("-pc_patch_construct_type", "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr);
1893   if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);}
1894   ierr = PetscOptionsInt("-pc_patch_vanka_dim", "Topological dimension of entities for Vanka to ignore", "PCPATCH", patch->vankadim, &patch->vankadim, &flg);CHKERRQ(ierr);
1895   ierr = PetscOptionsInt("-pc_patch_ignore_dim", "Topological dimension of entities for completion to ignore", "PCPATCH", patch->ignoredim, &patch->ignoredim, &flg);CHKERRQ(ierr);
1896   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);
1897   if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);}
1898   ierr = PetscOptionsBool("-pc_patch_symmetrise_sweep", "Go start->end, end->start?", "PCPATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr);
1899 
1900   /* If the user has set the number of subspaces, use that for the buffer size,
1901      otherwise use a large number */
1902   if (patch->nsubspaces <= 0) {
1903     nfields = 128;
1904   } else {
1905     nfields = patch->nsubspaces;
1906   }
1907   ierr = PetscMalloc1(nfields, &ifields);CHKERRQ(ierr);
1908   ierr = PetscOptionsGetIntArray(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_patch_exclude_subspaces",ifields,&nfields,&flg);CHKERRQ(ierr);
1909   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");
1910   if (flg) {
1911     PetscHSetIClear(patch->subspaces_to_exclude);
1912     for (k = 0; k < nfields; k++) {
1913       PetscHSetIAdd(patch->subspaces_to_exclude, ifields[k]);
1914     }
1915   }
1916   ierr = PetscFree(ifields);CHKERRQ(ierr);
1917 
1918   ierr = PetscOptionsBool("-pc_patch_patches_view", "Print out information during patch construction", "PCPATCH", patch->viewPatches, &patch->viewPatches, &flg);CHKERRQ(ierr);
1919   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, "-pc_patch_cells_view",   &patch->viewerCells,   &patch->formatCells,   &patch->viewCells);CHKERRQ(ierr);
1920   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, "-pc_patch_points_view",  &patch->viewerPoints,  &patch->formatPoints,  &patch->viewPoints);CHKERRQ(ierr);
1921   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, "-pc_patch_section_view", &patch->viewerSection, &patch->formatSection, &patch->viewSection);CHKERRQ(ierr);
1922   ierr = PetscOptionsGetViewer(comm,((PetscObject) pc)->options,prefix, "-pc_patch_mat_view",     &patch->viewerMatrix,  &patch->formatMatrix,  &patch->viewMatrix);CHKERRQ(ierr);
1923   ierr = PetscOptionsTail();CHKERRQ(ierr);
1924   patch->optionsSet = PETSC_TRUE;
1925   PetscFunctionReturn(0);
1926 }
1927 
1928 static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc)
1929 {
1930   PC_PATCH          *patch = (PC_PATCH*) pc->data;
1931   KSPConvergedReason reason;
1932   PetscInt           i;
1933   PetscErrorCode     ierr;
1934 
1935   PetscFunctionBegin;
1936   if (!patch->save_operators) {
1937     /* Can't do this here because the sub KSPs don't have an operator attached yet. */
1938     PetscFunctionReturn(0);
1939   }
1940   for (i = 0; i < patch->npatch; ++i) {
1941     if (!((KSP) patch->solver[i])->setfromoptionscalled) {
1942       ierr = KSPSetFromOptions((KSP) patch->solver[i]);CHKERRQ(ierr);
1943     }
1944     ierr = KSPSetUp((KSP) patch->solver[i]);CHKERRQ(ierr);
1945     ierr = KSPGetConvergedReason((KSP) patch->solver[i], &reason);CHKERRQ(ierr);
1946     if (reason == KSP_DIVERGED_PCSETUP_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1947   }
1948   PetscFunctionReturn(0);
1949 }
1950 
1951 static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer)
1952 {
1953   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1954   PetscViewer    sviewer;
1955   PetscBool      isascii;
1956   PetscMPIInt    rank;
1957   PetscErrorCode ierr;
1958 
1959   PetscFunctionBegin;
1960   /* TODO Redo tabbing with set tbas in new style */
1961   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
1962   if (!isascii) PetscFunctionReturn(0);
1963   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr);
1964   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1965   ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr);
1966   ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr);
1967   if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);}
1968   else                           {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);}
1969   if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);}
1970   else                         {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);}
1971   if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);}
1972   else                        {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);}
1973   if (patch->patchconstructop == PCPatchConstruct_Star)       {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);}
1974   else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);}
1975   else if (patch->patchconstructop == PCPatchConstruct_User)  {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);}
1976   else                                                        {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);}
1977   ierr = PetscViewerASCIIPrintf(viewer, "Solver on patches (all same):\n");CHKERRQ(ierr);
1978   if (patch->solver) {
1979     ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
1980     if (!rank) {
1981       ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr);
1982       ierr = PetscObjectView(patch->solver[0], sviewer);CHKERRQ(ierr);
1983       ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr);
1984     }
1985     ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
1986   } else {
1987     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1988     ierr = PetscViewerASCIIPrintf(viewer, "Solver not yet set.\n");CHKERRQ(ierr);
1989     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1990   }
1991   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1992   PetscFunctionReturn(0);
1993 }
1994 
1995 /*MC
1996   PCPATCH - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping
1997             small block additive preconditioners. Block definition is based on topology from
1998             a DM and equation numbering from a PetscSection.
1999 
2000   Options Database Keys:
2001 + -pc_patch_cells_view   - Views the process local cell numbers for each patch
2002 . -pc_patch_points_view  - Views the process local mesh point numbers for each patch
2003 . -pc_patch_g2l_view     - Views the map between global dofs and patch local dofs for each patch
2004 . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary
2005 - -pc_patch_sub_mat_view - Views the matrix associated with each patch
2006 
2007   Level: intermediate
2008 
2009 .seealso: PCType, PCCreate(), PCSetType()
2010 M*/
2011 PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc)
2012 {
2013   PC_PATCH      *patch;
2014   PetscErrorCode ierr;
2015 
2016   PetscFunctionBegin;
2017   ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr);
2018 
2019   if (patch->subspaces_to_exclude) {
2020     PetscHSetIDestroy(&patch->subspaces_to_exclude);
2021   }
2022   PetscHSetICreate(&patch->subspaces_to_exclude);
2023 
2024   /* Set some defaults */
2025   patch->combined           = PETSC_FALSE;
2026   patch->save_operators     = PETSC_TRUE;
2027   patch->partition_of_unity = PETSC_FALSE;
2028   patch->codim              = -1;
2029   patch->dim                = -1;
2030   patch->vankadim           = -1;
2031   patch->ignoredim          = -1;
2032   patch->patchconstructop   = PCPatchConstruct_Star;
2033   patch->symmetrise_sweep   = PETSC_FALSE;
2034   patch->npatch             = 0;
2035   patch->userIS             = NULL;
2036   patch->optionsSet         = PETSC_FALSE;
2037   patch->iterationSet       = NULL;
2038   patch->user_patches       = PETSC_FALSE;
2039   ierr = PetscStrallocpy(MATDENSE, (char **) &patch->sub_mat_type);CHKERRQ(ierr);
2040   patch->viewPatches        = PETSC_FALSE;
2041   patch->viewCells          = PETSC_FALSE;
2042   patch->viewPoints         = PETSC_FALSE;
2043   patch->viewSection        = PETSC_FALSE;
2044   patch->viewMatrix         = PETSC_FALSE;
2045   patch->setupsolver        = PCSetUp_PATCH_Linear;
2046   patch->applysolver        = PCApply_PATCH_Linear;
2047   patch->resetsolver        = PCReset_PATCH_Linear;
2048   patch->destroysolver      = PCDestroy_PATCH_Linear;
2049 
2050   pc->data                 = (void *) patch;
2051   pc->ops->apply           = PCApply_PATCH;
2052   pc->ops->applytranspose  = 0; /* PCApplyTranspose_PATCH; */
2053   pc->ops->setup           = PCSetUp_PATCH;
2054   pc->ops->reset           = PCReset_PATCH;
2055   pc->ops->destroy         = PCDestroy_PATCH;
2056   pc->ops->setfromoptions  = PCSetFromOptions_PATCH;
2057   pc->ops->setuponblocks   = PCSetUpOnBlocks_PATCH;
2058   pc->ops->view            = PCView_PATCH;
2059   pc->ops->applyrichardson = 0;
2060 
2061   PetscFunctionReturn(0);
2062 }
2063