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