xref: /petsc/src/ksp/pc/impls/patch/pcpatch.c (revision 5f824522500cf28ac329aa5a95292b3adbaf86e1)
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   PetscFunctionBeginHot;
14   ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
15   ierr = PetscObjectView(obj, viewer);CHKERRQ(ierr);
16   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
17   PetscFunctionReturn(0);
18 }
19 
20 static PetscErrorCode PCPatchConstruct_Star(void *vpatch, DM dm, PetscInt point, PetscHashI ht)
21 {
22   PetscInt       starSize;
23   PetscInt      *star = NULL, si;
24   PetscErrorCode ierr;
25 
26   PetscFunctionBegin;
27   PetscHashIClear(ht);
28   /* To start with, add the point we care about */
29   PetscHashIAdd(ht, point, 0);
30   /* Loop over all the points that this point connects to */
31   ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
32   for (si = 0; si < starSize*2; si += 2) {PetscHashIAdd(ht, star[si], 0);}
33   ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 static PetscErrorCode PCPatchConstruct_Vanka(void *vpatch, DM dm, PetscInt point, PetscHashI ht)
38 {
39   PC_PATCH      *patch = (PC_PATCH *) vpatch;
40   PetscInt       starSize;
41   PetscInt      *star = NULL;
42   PetscBool      shouldIgnore = PETSC_FALSE;
43   PetscInt       cStart, cEnd, iStart, iEnd, si;
44   PetscErrorCode ierr;
45 
46   PetscFunctionBegin;
47   PetscHashIClear(ht);
48   /* To start with, add the point we care about */
49   PetscHashIAdd(ht, point, 0);
50   /* Should we ignore any points of a certain dimension? */
51   if (patch->vankadim >= 0) {
52     shouldIgnore = PETSC_TRUE;
53     ierr = DMPlexGetDepthStratum(dm, patch->vankadim, &iStart, &iEnd);CHKERRQ(ierr);
54   }
55   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
56   /* Loop over all the cells that this point connects to */
57   ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
58   for (si = 0; si < starSize*2; si += 2) {
59     const PetscInt cell = star[si];
60     PetscInt       closureSize;
61     PetscInt      *closure = NULL, ci;
62 
63     if (cell < cStart || cell >= cEnd) continue;
64     /* now loop over all entities in the closure of that cell */
65     ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
66     for (ci = 0; ci < closureSize*2; ci += 2) {
67       const PetscInt newpoint = closure[ci];
68 
69       /* We've been told to ignore entities of this type.*/
70       if (shouldIgnore && newpoint >= iStart && newpoint < iEnd) continue;
71       PetscHashIAdd(ht, newpoint, 0);
72     }
73     ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
74   }
75   ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 /* The user's already set the patches in patch->userIS. Build the hash tables */
80 static PetscErrorCode PCPatchConstruct_User(void *vpatch, DM dm, PetscInt point, PetscHashI ht)
81 {
82   PC_PATCH       *patch   = (PC_PATCH *) vpatch;
83   IS              patchis = patch->userIS[point];
84   PetscInt        n;
85   const PetscInt *patchdata;
86   PetscInt        pStart, pEnd, i;
87   PetscErrorCode  ierr;
88 
89   PetscFunctionBegin;
90   PetscHashIClear(ht);
91   ierr = DMPlexGetChart(dm, &pStart, &pEnd);
92   ierr = ISGetLocalSize(patchis, &n);CHKERRQ(ierr);
93   ierr = ISGetIndices(patchis, &patchdata);CHKERRQ(ierr);
94   for (i = 0; i < n; ++i) {
95     const PetscInt ownedpoint = patchdata[i];
96 
97     if (ownedpoint < pStart || ownedpoint >= pEnd) {
98       SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D was not in [%D, %D)", ownedpoint, pStart, pEnd);
99     }
100     PetscHashIAdd(ht, ownedpoint, 0);
101   }
102   ierr = ISRestoreIndices(patchis, &patchdata);CHKERRQ(ierr);
103   PetscFunctionReturn(0);
104 }
105 
106 static PetscErrorCode PCPatchCreateDefaultSF_Private(PC pc, PetscInt n, const PetscSF *sf, const PetscInt *bs)
107 {
108   PC_PATCH      *patch = (PC_PATCH *) pc->data;
109   PetscInt       i;
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   if (n == 1 && bs[0] == 1) {
114     patch->defaultSF = sf[0];
115     ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr);
116   } else {
117     PetscInt     allRoots = 0, allLeaves = 0;
118     PetscInt     leafOffset = 0;
119     PetscInt    *ilocal = NULL;
120     PetscSFNode *iremote = NULL;
121     PetscInt    *remoteOffsets = NULL;
122     PetscInt     index = 0;
123     PetscHashI   rankToIndex;
124     PetscInt     numRanks = 0;
125     PetscSFNode *remote = NULL;
126     PetscSF      rankSF;
127     PetscInt    *ranks = NULL;
128     PetscInt    *offsets = NULL;
129     MPI_Datatype contig;
130     PetscHashI   ranksUniq;
131 
132     /* First figure out how many dofs there are in the concatenated numbering.
133      * allRoots: number of owned global dofs;
134      * allLeaves: number of visible dofs (global + ghosted).
135      */
136     for (i = 0; i < n; ++i) {
137       PetscInt nroots, nleaves;
138 
139       ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
140       allRoots  += nroots * bs[i];
141       allLeaves += nleaves * bs[i];
142     }
143     ierr = PetscMalloc1(allLeaves, &ilocal);CHKERRQ(ierr);
144     ierr = PetscMalloc1(allLeaves, &iremote);CHKERRQ(ierr);
145     /* Now build an SF that just contains process connectivity. */
146     PetscHashICreate(ranksUniq);
147     for (i = 0; i < n; ++i) {
148       const PetscMPIInt *ranks = NULL;
149       PetscInt           nranks, j;
150 
151       ierr = PetscSFSetUp(sf[i]);CHKERRQ(ierr);
152       ierr = PetscSFGetRanks(sf[i], &nranks, &ranks, NULL, NULL, NULL);CHKERRQ(ierr);
153       /* These are all the ranks who communicate with me. */
154       for (j = 0; j < nranks; ++j) {
155         PetscHashIAdd(ranksUniq, (PetscInt) ranks[j], 0);
156       }
157     }
158     PetscHashISize(ranksUniq, numRanks);
159     ierr = PetscMalloc1(numRanks, &remote);CHKERRQ(ierr);
160     ierr = PetscMalloc1(numRanks, &ranks);CHKERRQ(ierr);
161     PetscHashIGetKeys(ranksUniq, &index, ranks);
162 
163     PetscHashICreate(rankToIndex);
164     for (i = 0; i < numRanks; ++i) {
165       remote[i].rank  = ranks[i];
166       remote[i].index = 0;
167       PetscHashIAdd(rankToIndex, ranks[i], i);
168     }
169     ierr = PetscFree(ranks);CHKERRQ(ierr);
170     PetscHashIDestroy(ranksUniq);
171     ierr = PetscSFCreate(PetscObjectComm((PetscObject) pc), &rankSF);CHKERRQ(ierr);
172     ierr = PetscSFSetGraph(rankSF, 1, numRanks, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
173     ierr = PetscSFSetUp(rankSF);CHKERRQ(ierr);
174     /* OK, use it to communicate the root offset on the remote
175      * processes for each subspace. */
176     ierr = PetscMalloc1(n, &offsets);CHKERRQ(ierr);
177     ierr = PetscMalloc1(n*numRanks, &remoteOffsets);CHKERRQ(ierr);
178 
179     offsets[0] = 0;
180     for (i = 1; i < n; ++i) {
181       PetscInt nroots;
182 
183       ierr = PetscSFGetGraph(sf[i-1], &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
184       offsets[i] = offsets[i-1] + nroots*bs[i-1];
185     }
186     /* Offsets are the offsets on the current process of the
187      * global dof numbering for the subspaces. */
188     ierr = MPI_Type_contiguous(n, MPIU_INT, &contig);CHKERRQ(ierr);
189     ierr = MPI_Type_commit(&contig);CHKERRQ(ierr);
190 
191     ierr = PetscSFBcastBegin(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr);
192     ierr = PetscSFBcastEnd(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr);
193     ierr = MPI_Type_free(&contig);CHKERRQ(ierr);
194     ierr = PetscFree(offsets);CHKERRQ(ierr);
195     ierr = PetscSFDestroy(&rankSF);CHKERRQ(ierr);
196     /* Now remoteOffsets contains the offsets on the remote
197      * processes who communicate with me.  So now we can
198      * concatenate the list of SFs into a single one. */
199     index = 0;
200     for (i = 0; i < n; ++i) {
201       const PetscSFNode *remote = NULL;
202       const PetscInt    *local  = NULL;
203       PetscInt           nroots, nleaves, j;
204 
205       ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, &local, &remote);CHKERRQ(ierr);
206       for (j = 0; j < nleaves; ++j) {
207         PetscInt rank = remote[j].rank;
208         PetscInt idx, rootOffset, k;
209 
210         PetscHashIMap(rankToIndex, rank, idx);
211         if (idx == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Didn't find rank, huh?");
212         /* Offset on given rank for ith subspace */
213         rootOffset = remoteOffsets[n*idx + i];
214         for (k = 0; k < bs[i]; ++k) {
215           ilocal[index]        = local[j]*bs[i] + k + leafOffset;
216           iremote[index].rank  = remote[j].rank;
217           iremote[index].index = remote[j].index*bs[i] + k + rootOffset;
218           ++index;
219         }
220       }
221       leafOffset += nleaves * bs[i];
222     }
223     PetscHashIDestroy(rankToIndex);
224     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
225     ierr = PetscSFCreate(PetscObjectComm((PetscObject)pc), &patch->defaultSF);CHKERRQ(ierr);
226     ierr = PetscSFSetGraph(patch->defaultSF, allRoots, allLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
227   }
228   PetscFunctionReturn(0);
229 }
230 
231 /* TODO: Docs */
232 PetscErrorCode PCPatchSetIgnoreDim(PC pc, PetscInt dim)
233 {
234   PC_PATCH *patch = (PC_PATCH *) pc->data;
235   PetscFunctionBegin;
236   patch->ignoredim = dim;
237   PetscFunctionReturn(0);
238 }
239 
240 /* TODO: Docs */
241 PetscErrorCode PCPatchGetIgnoreDim(PC pc, PetscInt *dim)
242 {
243   PC_PATCH *patch = (PC_PATCH *) pc->data;
244   PetscFunctionBegin;
245   *dim = patch->ignoredim;
246   PetscFunctionReturn(0);
247 }
248 
249 /* TODO: Docs */
250 PetscErrorCode PCPatchSetSaveOperators(PC pc, PetscBool flg)
251 {
252   PC_PATCH *patch = (PC_PATCH *) pc->data;
253   PetscFunctionBegin;
254   patch->save_operators = flg;
255   PetscFunctionReturn(0);
256 }
257 
258 /* TODO: Docs */
259 PetscErrorCode PCPatchGetSaveOperators(PC pc, PetscBool *flg)
260 {
261   PC_PATCH *patch = (PC_PATCH *) pc->data;
262   PetscFunctionBegin;
263   *flg = patch->save_operators;
264   PetscFunctionReturn(0);
265 }
266 
267 /* TODO: Docs */
268 PetscErrorCode PCPatchSetPartitionOfUnity(PC pc, PetscBool flg)
269 {
270     PC_PATCH *patch = (PC_PATCH *) pc->data;
271     PetscFunctionBegin;
272     patch->partition_of_unity = flg;
273     PetscFunctionReturn(0);
274 }
275 
276 /* TODO: Docs */
277 PetscErrorCode PCPatchGetPartitionOfUnity(PC pc, PetscBool *flg)
278 {
279     PC_PATCH *patch = (PC_PATCH *) pc->data;
280     PetscFunctionBegin;
281     *flg = patch->partition_of_unity;
282     PetscFunctionReturn(0);
283 }
284 
285 /* TODO: Docs */
286 PetscErrorCode PCPatchSetMultiplicative(PC pc, PetscBool flg)
287 {
288   PC_PATCH *patch = (PC_PATCH *) pc->data;
289   PetscFunctionBegin;
290   patch->multiplicative = flg;
291   PetscFunctionReturn(0);
292 }
293 
294 /* TODO: Docs */
295 PetscErrorCode PCPatchGetMultiplicative(PC pc, PetscBool *flg)
296 {
297   PC_PATCH *patch = (PC_PATCH *) pc->data;
298   PetscFunctionBegin;
299   *flg = patch->multiplicative;
300   PetscFunctionReturn(0);
301 }
302 
303 /* TODO: Docs */
304 PetscErrorCode PCPatchSetSubMatType(PC pc, MatType sub_mat_type)
305 {
306   PC_PATCH      *patch = (PC_PATCH *) pc->data;
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   if (patch->sub_mat_type) {ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);}
311   ierr = PetscStrallocpy(sub_mat_type, (char **) &patch->sub_mat_type);CHKERRQ(ierr);
312   PetscFunctionReturn(0);
313 }
314 
315 /* TODO: Docs */
316 PetscErrorCode PCPatchGetSubMatType(PC pc, MatType *sub_mat_type)
317 {
318   PC_PATCH *patch = (PC_PATCH *) pc->data;
319   PetscFunctionBegin;
320   *sub_mat_type = patch->sub_mat_type;
321   PetscFunctionReturn(0);
322 }
323 
324 /* TODO: Docs */
325 PetscErrorCode PCPatchSetCellNumbering(PC pc, PetscSection cellNumbering)
326 {
327   PC_PATCH      *patch = (PC_PATCH *) pc->data;
328   PetscErrorCode ierr;
329 
330   PetscFunctionBegin;
331   patch->cellNumbering = cellNumbering;
332   ierr = PetscObjectReference((PetscObject) cellNumbering);CHKERRQ(ierr);
333   PetscFunctionReturn(0);
334 }
335 
336 /* TODO: Docs */
337 PetscErrorCode PCPatchGetCellNumbering(PC pc, PetscSection *cellNumbering)
338 {
339   PC_PATCH *patch = (PC_PATCH *) pc->data;
340   PetscFunctionBegin;
341   *cellNumbering = patch->cellNumbering;
342   PetscFunctionReturn(0);
343 }
344 
345 /* TODO: Docs */
346 PetscErrorCode PCPatchSetConstructType(PC pc, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx)
347 {
348   PC_PATCH *patch = (PC_PATCH *) pc->data;
349 
350   PetscFunctionBegin;
351   patch->ctype = ctype;
352   switch (ctype) {
353   case PC_PATCH_STAR:
354     patch->patchconstructop = PCPatchConstruct_Star;
355     break;
356   case PC_PATCH_VANKA:
357     patch->patchconstructop = PCPatchConstruct_Vanka;
358     break;
359   case PC_PATCH_USER:
360   case PC_PATCH_PYTHON:
361     patch->user_patches            = PETSC_TRUE;
362     patch->patchconstructop        = PCPatchConstruct_User;
363     patch->userpatchconstructionop = func;
364     patch->userpatchconstructctx   = ctx;
365     break;
366   default:
367     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype);
368   }
369   PetscFunctionReturn(0);
370 }
371 
372 /* TODO: Docs */
373 PetscErrorCode PCPatchGetConstructType(PC pc, PCPatchConstructType *ctype, PetscErrorCode (**func)(PC, PetscInt *, IS **, IS *, void *), void **ctx)
374 {
375   PC_PATCH *patch = (PC_PATCH *) pc->data;
376 
377   PetscFunctionBegin;
378   *ctype = patch->ctype;
379   switch (patch->ctype) {
380   case PC_PATCH_STAR:
381   case PC_PATCH_VANKA:
382     break;
383   case PC_PATCH_USER:
384   case PC_PATCH_PYTHON:
385     *func = patch->userpatchconstructionop;
386     *ctx  = patch->userpatchconstructctx;
387     break;
388   default:
389     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype);
390   }
391   PetscFunctionReturn(0);
392 }
393 
394 /* TODO: Docs */
395 PetscErrorCode PCPatchSetDiscretisationInfo(PC pc, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap,
396                                             const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes)
397 {
398   PC_PATCH      *patch = (PC_PATCH *) pc->data;
399   DM             dm;
400   PetscSF       *sfs;
401   PetscInt       cStart, cEnd, i, j;
402   PetscErrorCode ierr;
403 
404   PetscFunctionBegin;
405   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
406   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
407   ierr = PetscMalloc1(nsubspaces, &sfs);CHKERRQ(ierr);
408   ierr = PetscMalloc1(nsubspaces, &patch->dofSection);CHKERRQ(ierr);
409   ierr = PetscMalloc1(nsubspaces, &patch->bs);CHKERRQ(ierr);
410   ierr = PetscMalloc1(nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr);
411   ierr = PetscMalloc1(nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr);
412   ierr = PetscMalloc1(nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr);
413 
414   patch->nsubspaces       = nsubspaces;
415   patch->totalDofsPerCell = 0;
416   for (i = 0; i < nsubspaces; ++i) {
417     ierr = DMGetDefaultSection(dms[i], &patch->dofSection[i]);CHKERRQ(ierr);
418     ierr = PetscObjectReference((PetscObject) patch->dofSection[i]);CHKERRQ(ierr);
419     ierr = DMGetDefaultSF(dms[i], &sfs[i]);CHKERRQ(ierr);
420     patch->bs[i]              = bs[i];
421     patch->nodesPerCell[i]    = nodesPerCell[i];
422     patch->totalDofsPerCell  += nodesPerCell[i]*bs[i];
423     ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i]*bs[i], &patch->cellNodeMap[i]);CHKERRQ(ierr);
424     for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]*bs[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j];
425     patch->subspaceOffsets[i] = subspaceOffsets[i];
426   }
427   ierr = PCPatchCreateDefaultSF_Private(pc, nsubspaces, sfs, patch->bs);CHKERRQ(ierr);
428   ierr = PetscFree(sfs);CHKERRQ(ierr);
429 
430   patch->subspaceOffsets[nsubspaces] = subspaceOffsets[nsubspaces];
431   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr);
432   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr);
433   PetscFunctionReturn(0);
434 }
435 
436 /* TODO: Docs */
437 PetscErrorCode PCPatchSetDiscretisationInfoCombined(PC pc, DM dm, PetscInt *nodesPerCell, const PetscInt **cellNodeMap, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes)
438 {
439   PC_PATCH      *patch = (PC_PATCH *) pc->data;
440   PetscInt       cStart, cEnd, i, j;
441   PetscErrorCode ierr;
442 
443   PetscFunctionBegin;
444   patch->combined = PETSC_TRUE;
445   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
446   ierr = DMGetNumFields(dm, &patch->nsubspaces);CHKERRQ(ierr);
447   ierr = PetscCalloc1(patch->nsubspaces, &patch->dofSection);CHKERRQ(ierr);
448   ierr = PetscMalloc1(patch->nsubspaces, &patch->bs);CHKERRQ(ierr);
449   ierr = PetscMalloc1(patch->nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr);
450   ierr = PetscMalloc1(patch->nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr);
451   ierr = PetscCalloc1(patch->nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr);
452   ierr = DMGetDefaultSection(dm, &patch->dofSection[0]);CHKERRQ(ierr);
453   ierr = PetscObjectReference((PetscObject) patch->dofSection[0]);CHKERRQ(ierr);
454   ierr = PetscSectionGetStorageSize(patch->dofSection[0], &patch->subspaceOffsets[patch->nsubspaces]);CHKERRQ(ierr);
455   patch->totalDofsPerCell = 0;
456   for (i = 0; i < patch->nsubspaces; ++i) {
457     patch->bs[i]             = 1;
458     patch->nodesPerCell[i]   = nodesPerCell[i];
459     patch->totalDofsPerCell += nodesPerCell[i];
460     ierr = PetscMalloc1((cEnd-cStart)*nodesPerCell[i], &patch->cellNodeMap[i]);CHKERRQ(ierr);
461     for (j = 0; j < (cEnd-cStart)*nodesPerCell[i]; ++j) patch->cellNodeMap[i][j] = cellNodeMap[i][j];
462   }
463   ierr = DMGetDefaultSF(dm, &patch->defaultSF);CHKERRQ(ierr);
464   ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr);
465   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr);
466   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr);
467   PetscFunctionReturn(0);
468 }
469 
470 /*@C
471 
472   PCPatchSetComputeOperator - Set the callback used to compute patch matrices
473 
474   Input Parameters:
475 + pc   - The PC
476 . func - The callback
477 - ctx  - The user context
478 
479   Level: advanced
480 
481   Note:
482   The callback has signature:
483 +  usercomputeop(pc, mat, ncell, cells, n, u, ctx)
484 +  pc    - The PC
485 +  mat   - The patch matrix
486 +  ncell - The number of cells to integrate over
487 +  cells - An array of the cell numbers
488 +  n     - The size of g2l
489 +  g2l   - The global to local dof translation table
490 +  ctx   - The user context
491   and can assume that the matrix entries have been set to zero before the call.
492 
493 .seealso: PCPatchGetComputeOperator(), PCPatchSetDiscretisationInfo()
494 @*/
495 PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, PetscInt, Mat, PetscInt, const PetscInt [], PetscInt, const PetscInt *, void *), void *ctx)
496 {
497   PC_PATCH *patch = (PC_PATCH *) pc->data;
498 
499   PetscFunctionBegin;
500   patch->usercomputeop  = func;
501   patch->usercomputectx = ctx;
502   PetscFunctionReturn(0);
503 }
504 
505 /* On entry, ht contains the topological entities whose dofs we are responsible for solving for;
506    on exit, cht contains all the topological entities we need to compute their residuals.
507    In full generality this should incorporate knowledge of the sparsity pattern of the matrix;
508    here we assume a standard FE sparsity pattern.*/
509 /* TODO: Use DMPlexGetAdjacency() */
510 /* TODO: Look at temp buffer management for GetClosure() */
511 static PetscErrorCode PCPatchCompleteCellPatch(PC pc, PetscHashI ht, PetscHashI cht)
512 {
513   DM             dm;
514   PetscHashIIter hi;
515   PetscInt       point;
516   PetscInt      *star = NULL, *closure = NULL;
517   PetscInt       ignoredim, iStart, iEnd, starSize, closureSize, si, ci;
518   PetscErrorCode ierr;
519 
520   PetscFunctionBegin;
521   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
522   ierr = PCPatchGetIgnoreDim(pc, &ignoredim);CHKERRQ(ierr);
523   if (ignoredim >= 0) {ierr = DMPlexGetDepthStratum(dm, ignoredim, &iStart, &iEnd);CHKERRQ(ierr);}
524   PetscHashIClear(cht);
525   PetscHashIIterBegin(ht, hi);
526   while (!PetscHashIIterAtEnd(ht, hi)) {
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       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         PetscInt seenpoint = closure[ci];
539         if (ignoredim >= 0 && seenpoint >= iStart && seenpoint < iEnd) continue;
540         PetscHashIAdd(cht, seenpoint, 0);
541       }
542     }
543     ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr);
544   }
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   ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr);
1264   ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr);
1265   if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);}
1266   /* Sparse patch matrices */
1267   if (!flg) {
1268     PetscBT         bt;
1269     PetscInt       *dnnz      = NULL;
1270     const PetscInt *dofsArray = NULL;
1271     PetscInt        pStart, pEnd, ncell, offset, c, i, j;
1272 
1273     ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1274     ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1275     point += pStart;
1276     if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1277     ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1278     ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1279     ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr);
1280     ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
1281     /* XXX: This uses N^2 bits to store the sparsity pattern on a
1282      * patch.  This is probably OK if the patches are not too big,
1283      * but could use quite a bit of memory for planes in 3D.
1284      * Should we switch based on the value of rsize to a
1285      * hash-table (slower, but more memory efficient) approach? */
1286     ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr);
1287     for (c = 0; c < ncell; ++c) {
1288       const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell;
1289       for (i = 0; i < patch->totalDofsPerCell; ++i) {
1290         const PetscInt row = idx[i];
1291         for (j = 0; j < patch->totalDofsPerCell; ++j) {
1292           const PetscInt col = idx[j];
1293           const PetscInt key = row*rsize + col;
1294           if (!PetscBTLookupSet(bt, key)) ++dnnz[row];
1295         }
1296       }
1297     }
1298     ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
1299     ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr);
1300     ierr = PetscFree(dnnz);CHKERRQ(ierr);
1301     ierr = PCPatchZeroFillMatrix_Private(*mat, ncell, patch->totalDofsPerCell, &dofsArray[offset*patch->totalDofsPerCell]);CHKERRQ(ierr);
1302     ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
1303     ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1304   }
1305   ierr = MatSetUp(*mat);CHKERRQ(ierr);
1306   PetscFunctionReturn(0);
1307 }
1308 
1309 static PetscErrorCode PCPatchComputeOperator_DMPlex_Private(PC pc, PetscInt patchNum, Mat J, PetscInt ncell, const PetscInt cells[], PetscInt n, const PetscInt *l2p, void *ctx)
1310 {
1311   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1312   DM              dm;
1313   PetscSection    s;
1314   const PetscInt *parray, *oarray;
1315   PetscInt        Nf = patch->nsubspaces, Np, poff, p, f;
1316   PetscErrorCode  ierr;
1317 
1318   PetscFunctionBegin;
1319   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1320   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1321   /* Set offset into patch */
1322   ierr = PetscSectionGetDof(patch->pointCounts, patchNum, &Np);CHKERRQ(ierr);
1323   ierr = PetscSectionGetOffset(patch->pointCounts, patchNum, &poff);CHKERRQ(ierr);
1324   ierr = ISGetIndices(patch->points, &parray);CHKERRQ(ierr);
1325   ierr = ISGetIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1326   for (f = 0; f < Nf; ++f) {
1327     for (p = 0; p < Np; ++p) {
1328       const PetscInt point = parray[poff+p];
1329       PetscInt       dof;
1330 
1331       ierr = PetscSectionGetFieldDof(patch->patchSection, point, f, &dof);CHKERRQ(ierr);
1332       ierr = PetscSectionSetFieldOffset(patch->patchSection, point, f, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);
1333       if (patch->nsubspaces == 1) {ierr = PetscSectionSetOffset(patch->patchSection, point, oarray[(poff+p)*Nf+f]);CHKERRQ(ierr);}
1334       else                        {ierr = PetscSectionSetOffset(patch->patchSection, point, -1);CHKERRQ(ierr);}
1335     }
1336   }
1337   ierr = ISRestoreIndices(patch->points, &parray);CHKERRQ(ierr);
1338   ierr = ISRestoreIndices(patch->offs,   &oarray);CHKERRQ(ierr);
1339   if (patch->viewSection) {ierr = ObjectView((PetscObject) patch->patchSection, patch->viewerSection, patch->formatSection);CHKERRQ(ierr);}
1340   /* TODO Shut off MatViewFromOptions() in MatAssemblyEnd() here */
1341   ierr = DMPlexComputeJacobian_Patch_Internal(pc->dm, patch->patchSection, patch->patchSection, 0, ncell, cells, 0.0, 0.0, NULL, NULL, J, J, ctx);CHKERRQ(ierr);
1342   PetscFunctionReturn(0);
1343 }
1344 
1345 static PetscErrorCode PCPatchComputeOperator_Private(PC pc, Mat mat, Mat multMat, PetscInt point)
1346 {
1347   PC_PATCH       *patch = (PC_PATCH *) pc->data;
1348   const PetscInt *dofsArray;
1349   const PetscInt *cellsArray;
1350   PetscInt        ncell, offset, pStart, pEnd;
1351   PetscErrorCode  ierr;
1352 
1353   PetscFunctionBegin;
1354   ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1355   if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n");
1356   ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1357   ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1358   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
1359 
1360   point += pStart;
1361   if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
1362 
1363   ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
1364   ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
1365   if (ncell <= 0) {
1366     ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1367     PetscFunctionReturn(0);
1368   }
1369   PetscStackPush("PCPatch user callback");
1370   ierr = patch->usercomputeop(pc, point, mat, ncell, cellsArray + offset, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputectx);CHKERRQ(ierr);
1371   PetscStackPop;
1372   ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
1373   ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
1374   /* TODO Can apply these BCs through the patchSection (would need to changes sizes of patchX/Y etc) */
1375   /* Apply boundary conditions.  Could also do this through the local_to_patch guy. */
1376   if (patch->multiplicative) {
1377     ierr = MatCopy(mat, multMat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1378     ierr = MatZeroRowsColumnsIS(multMat, patch->multBcs[point-pStart], (PetscScalar) 1.0, NULL, NULL);CHKERRQ(ierr);
1379   }
1380   ierr = MatZeroRowsColumnsIS(mat, patch->bcs[point-pStart], (PetscScalar) 1.0, NULL, NULL);CHKERRQ(ierr);
1381   if (patch->viewMatrix) {ierr = ObjectView((PetscObject) mat, patch->viewerMatrix, patch->formatMatrix);CHKERRQ(ierr);}
1382   ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
1383   PetscFunctionReturn(0);
1384 }
1385 
1386 static PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat)
1387 {
1388   PC_PATCH          *patch     = (PC_PATCH *) pc->data;
1389   const PetscScalar *xArray    = NULL;
1390   PetscScalar       *yArray    = NULL;
1391   const PetscInt    *gtolArray = NULL;
1392   PetscInt           dof, offset, lidx;
1393   PetscErrorCode     ierr;
1394 
1395   PetscFunctionBeginHot;
1396   ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
1397   ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr);
1398   ierr = VecGetArray(y, &yArray);CHKERRQ(ierr);
1399   ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
1400   ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr);
1401   ierr = ISGetIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1402   if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n");
1403   if (mode == ADD_VALUES    && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n");
1404   for (lidx = 0; lidx < dof; ++lidx) {
1405     const PetscInt gidx = gtolArray[offset+lidx];
1406 
1407     if (mode == INSERT_VALUES) yArray[lidx]  = xArray[gidx]; /* Forward */
1408     else                       yArray[gidx] += xArray[lidx]; /* Reverse */
1409   }
1410   ierr = ISRestoreIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
1411   ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr);
1412   ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr);
1413   ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
1414   PetscFunctionReturn(0);
1415 }
1416 
1417 static PetscErrorCode PCSetUp_PATCH(PC pc)
1418 {
1419   PC_PATCH       *patch   = (PC_PATCH *) pc->data;
1420   PetscScalar    *patchX  = NULL;
1421   const PetscInt *bcNodes = NULL;
1422   PetscInt        numBcs, i, j;
1423   const char     *prefix;
1424   PetscErrorCode  ierr;
1425 
1426   PetscFunctionBegin;
1427   if (!pc->setupcalled) {
1428     PetscInt pStart, pEnd, p;
1429     PetscInt localSize;
1430 
1431     ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
1432 
1433     if (!patch->nsubspaces) {
1434       DM           dm;
1435       PetscDS      prob;
1436       PetscSection s;
1437       PetscInt     cStart, cEnd, c, Nf, f, fsize, numGlobalBcs = 0, *globalBcs, *Nb, totNb = 0, **cellDofs;
1438 
1439       ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
1440       if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Must set DM for PCPATCH or call PCPatchSetDiscretisationInfo()");
1441       ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
1442       ierr = PetscSectionGetNumFields(s, &Nf);CHKERRQ(ierr);
1443       ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
1444       for (p = pStart; p < pEnd; ++p) {
1445         PetscInt cdof;
1446         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1447         numGlobalBcs += cdof;
1448       }
1449       ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1450       ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1451       ierr = PetscMalloc3(Nf, &Nb, Nf, &cellDofs, numGlobalBcs, &globalBcs);CHKERRQ(ierr);
1452       for (f = 0; f < Nf; ++f) {
1453         PetscFE        fe;
1454         PetscDualSpace sp;
1455         PetscSection   fs;
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 = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr);
1465         ierr = PetscSectionGetStorageSize(fs, &fsize);CHKERRQ(ierr);
1466         ierr = PetscMalloc1((cEnd-cStart)*Nb[f], &cellDofs[f]);CHKERRQ(ierr);
1467         for (c = cStart; c < cEnd; ++c) {
1468           PetscInt *closure = NULL;
1469           PetscInt  clSize  = 0, cl;
1470 
1471           ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
1472           for (cl = 0; cl < clSize*2; cl += 2) {
1473             const PetscInt p = closure[cl];
1474             PetscInt       fdof, d, foff;
1475 
1476             ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
1477             ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
1478             for (d = 0; d < fdof; ++d, ++cdoff) cellDofs[f][cdoff] = foff + d;
1479           }
1480           ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
1481         }
1482         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]);
1483       }
1484       numGlobalBcs = 0;
1485       for (p = pStart; p < pEnd; ++p) {
1486         const PetscInt *ind;
1487         PetscInt        off, cdof, d;
1488 
1489         ierr = PetscSectionGetOffset(s, p, &off);CHKERRQ(ierr);
1490         ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr);
1491         ierr = PetscSectionGetConstraintIndices(s, p, &ind);CHKERRQ(ierr);
1492         for (d = 0; d < cdof; ++d) globalBcs[numGlobalBcs++] = off + ind[d];
1493       }
1494 
1495       ierr = PCPatchSetDiscretisationInfoCombined(pc, dm, Nb, (const PetscInt **) cellDofs, numGlobalBcs, globalBcs, numGlobalBcs, globalBcs);CHKERRQ(ierr);
1496       for (f = 0; f < Nf; ++f) {
1497         ierr = PetscFree(cellDofs[f]);CHKERRQ(ierr);
1498       }
1499       ierr = PetscFree3(Nb, cellDofs, globalBcs);CHKERRQ(ierr);
1500       ierr = PCPatchSetComputeOperator(pc, PCPatchComputeOperator_DMPlex_Private, NULL);CHKERRQ(ierr);
1501     }
1502 
1503     localSize = patch->subspaceOffsets[patch->nsubspaces];
1504     ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localX);CHKERRQ(ierr);
1505     ierr = VecSetUp(patch->localX);CHKERRQ(ierr);
1506     ierr = VecDuplicate(patch->localX, &patch->localY);CHKERRQ(ierr);
1507     ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr);
1508     ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr);
1509     ierr = PCPatchCreateCellPatchBCs(pc);CHKERRQ(ierr);
1510 
1511     /* OK, now build the work vectors */
1512     ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr);
1513     ierr = PetscMalloc1(patch->npatch, &patch->patchX);CHKERRQ(ierr);
1514     ierr = PetscMalloc1(patch->npatch, &patch->patchY);CHKERRQ(ierr);
1515     if (patch->partition_of_unity && patch->multiplicative) {ierr = PetscMalloc1(patch->npatch, &patch->patch_dof_weights);CHKERRQ(ierr);}
1516     for (p = pStart; p < pEnd; ++p) {
1517       PetscInt dof;
1518 
1519       ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
1520       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchX[p-pStart]);CHKERRQ(ierr);
1521       ierr = VecSetUp(patch->patchX[p-pStart]);CHKERRQ(ierr);
1522       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchY[p-pStart]);CHKERRQ(ierr);
1523       ierr = VecSetUp(patch->patchY[p-pStart]);CHKERRQ(ierr);
1524       if (patch->partition_of_unity && patch->multiplicative) {
1525         ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patch_dof_weights[p-pStart]);CHKERRQ(ierr);
1526         ierr = VecSetUp(patch->patch_dof_weights[p-pStart]);CHKERRQ(ierr);
1527       }
1528     }
1529     ierr = PetscMalloc1(patch->npatch, &patch->ksp);CHKERRQ(ierr);
1530     ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
1531     for (i = 0; i < patch->npatch; ++i) {
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     }
1536     if (patch->save_operators) {
1537       ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr);
1538       if (patch->multiplicative) {ierr = PetscMalloc1(patch->npatch, &patch->multMat);CHKERRQ(ierr);}
1539       for (i = 0; i < patch->npatch; ++i) {
1540         ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i]);CHKERRQ(ierr);
1541         if (patch->multiplicative) {ierr = MatDuplicate(patch->mat[i], MAT_SHARE_NONZERO_PATTERN, &patch->multMat[i]);CHKERRQ(ierr);}
1542       }
1543     }
1544     ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
1545 
1546     /* If desired, calculate weights for dof multiplicity */
1547     if (patch->partition_of_unity) {
1548       ierr = VecDuplicate(patch->localX, &patch->dof_weights);CHKERRQ(ierr);
1549       for (i = 0; i < patch->npatch; ++i) {
1550         PetscInt dof;
1551 
1552         ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr);
1553         if (dof <= 0) continue;
1554         ierr = VecSet(patch->patchX[i], 1.0);CHKERRQ(ierr);
1555         /* TODO: Do we need different scatters for X and Y? */
1556         ierr = VecGetArray(patch->patchX[i], &patchX);CHKERRQ(ierr);
1557         /* Apply bcs to patchX (zero entries) */
1558         ierr = ISGetLocalSize(patch->bcs[i], &numBcs);CHKERRQ(ierr);
1559         ierr = ISGetIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr);
1560         for (j = 0; j < numBcs; ++j) patchX[bcNodes[j]] = 0;
1561         ierr = ISRestoreIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr);
1562         ierr = VecRestoreArray(patch->patchX[i], &patchX);CHKERRQ(ierr);
1563 
1564         ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchX[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
1565       }
1566       ierr = VecReciprocal(patch->dof_weights);CHKERRQ(ierr);
1567       if (patch->partition_of_unity && patch->multiplicative) {
1568         for (i = 0; i < patch->npatch; ++i) {
1569           ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->dof_weights, patch->patch_dof_weights[i], INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
1570         }
1571       }
1572     }
1573   }
1574   if (patch->save_operators) {
1575     for (i = 0; i < patch->npatch; ++i) {
1576       ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr);
1577       if (patch->multiplicative) {ierr = PCPatchComputeOperator_Private(pc, patch->mat[i], patch->multMat[i], i);CHKERRQ(ierr);}
1578       else                       {ierr = PCPatchComputeOperator_Private(pc, patch->mat[i], NULL, i);CHKERRQ(ierr);}
1579       ierr = KSPSetOperators(patch->ksp[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr);
1580     }
1581   }
1582   if (!pc->setupcalled && patch->optionsSet) for (i = 0; i < patch->npatch; ++i) {ierr = KSPSetFromOptions(patch->ksp[i]);CHKERRQ(ierr);}
1583   PetscFunctionReturn(0);
1584 }
1585 
1586 static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y)
1587 {
1588   PC_PATCH          *patch    = (PC_PATCH *) pc->data;
1589   const PetscScalar *globalX  = NULL;
1590   PetscScalar       *localX   = NULL;
1591   PetscScalar       *globalY  = NULL;
1592   PetscScalar       *patchX   = NULL;
1593   const PetscInt    *bcNodes  = NULL;
1594   PetscInt           nsweep   = patch->symmetrise_sweep ? 2 : 1;
1595   PetscInt           start[2] = {0, 0};
1596   PetscInt           end[2]   = {-1, -1};
1597   const PetscInt     inc[2]   = {1, -1};
1598   const PetscScalar *localY;
1599   const PetscInt    *iterationSet;
1600   PetscInt           pStart, numBcs, n, sweep, bc, j;
1601   PetscErrorCode     ierr;
1602 
1603   PetscFunctionBegin;
1604   ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
1605   ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr);
1606   end[0]   = patch->npatch;
1607   start[1] = patch->npatch-1;
1608   if (patch->user_patches) {
1609     ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr);
1610     start[1] = end[0] - 1;
1611     ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);
1612   }
1613   /* Scatter from global space into overlapped local spaces */
1614   ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr);
1615   ierr = VecGetArray(patch->localX, &localX);CHKERRQ(ierr);
1616   ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr);
1617   ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr);
1618   ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr);
1619   ierr = VecRestoreArray(patch->localX, &localX);CHKERRQ(ierr);
1620 
1621   ierr = VecSet(patch->localY, 0.0);CHKERRQ(ierr);
1622   ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr);
1623   for (sweep = 0; sweep < nsweep; sweep++) {
1624     for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) {
1625       Mat      multMat = NULL;
1626       PetscInt i       = patch->user_patches ? iterationSet[j] : j;
1627       PetscInt start, len;
1628 
1629       ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr);
1630       ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr);
1631       /* TODO: Squash out these guys in the setup as well. */
1632       if (len <= 0) continue;
1633       /* TODO: Do we need different scatters for X and Y? */
1634       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localX, patch->patchX[i], INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
1635       /* Apply bcs to patchX (zero entries) */
1636       ierr = VecGetArray(patch->patchX[i], &patchX);CHKERRQ(ierr);
1637       ierr = ISGetLocalSize(patch->bcs[i], &numBcs);CHKERRQ(ierr);
1638       ierr = ISGetIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr);
1639       for (bc = 0; bc < numBcs; ++bc) patchX[bcNodes[bc]] = 0;
1640       ierr = ISRestoreIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr);
1641       ierr = VecRestoreArray(patch->patchX[i], &patchX);CHKERRQ(ierr);
1642       if (!patch->save_operators) {
1643         Mat mat;
1644 
1645         ierr = PCPatchCreateMatrix_Private(pc, i, &mat);CHKERRQ(ierr);
1646         if (patch->multiplicative) {ierr = MatDuplicate(mat, MAT_SHARE_NONZERO_PATTERN, &multMat);CHKERRQ(ierr);}
1647         /* Populate operator here. */
1648         ierr = PCPatchComputeOperator_Private(pc, mat, multMat, i);CHKERRQ(ierr);
1649         ierr = KSPSetOperators(patch->ksp[i], mat, mat);
1650         /* Drop reference so the KSPSetOperators below will blow it away. */
1651         ierr = MatDestroy(&mat);CHKERRQ(ierr);
1652       } else if (patch->multiplicative) {
1653         multMat = patch->multMat[i];
1654       }
1655 
1656       ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1657       ierr = KSPSolve(patch->ksp[i], patch->patchX[i], patch->patchY[i]);CHKERRQ(ierr);
1658       ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
1659 
1660       if (!patch->save_operators) {
1661         PC pc;
1662         ierr = KSPSetOperators(patch->ksp[i], NULL, NULL);CHKERRQ(ierr);
1663         ierr = KSPGetPC(patch->ksp[i], &pc);CHKERRQ(ierr);
1664         /* Destroy PC context too, otherwise the factored matrix hangs around. */
1665         ierr = PCReset(pc);CHKERRQ(ierr);
1666       }
1667 
1668       if (patch->partition_of_unity && patch->multiplicative) {
1669         ierr = VecPointwiseMult(patch->patchY[i], patch->patchY[i], patch->patch_dof_weights[i]);CHKERRQ(ierr);
1670       }
1671       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchY[i], patch->localY, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
1672 
1673       /* Get the matrix on the patch but with only global bcs applied.
1674        * This matrix is then multiplied with the result from the previous solve
1675        * to obtain by how much the residual changes. */
1676       if (patch->multiplicative) {
1677         ierr = MatMult(multMat, patch->patchY[i], patch->patchX[i]);CHKERRQ(ierr);
1678         ierr = VecScale(patch->patchX[i], -1.0);CHKERRQ(ierr);
1679         ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchX[i], patch->localX, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
1680         if (!patch->save_operators) {ierr = MatDestroy(&multMat);CHKERRQ(ierr);}
1681       }
1682     }
1683   }
1684   if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);}
1685   /* XXX: should we do this on the global vector? */
1686   if (patch->partition_of_unity && !patch->multiplicative) {
1687     ierr = VecPointwiseMult(patch->localY, patch->localY, patch->dof_weights);CHKERRQ(ierr);
1688   }
1689   /* Now patch->localY contains the solution of the patch solves, so we need to combine them all. */
1690   ierr = VecSet(y, 0.0);CHKERRQ(ierr);
1691   ierr = VecGetArray(y, &globalY);CHKERRQ(ierr);
1692   ierr = VecGetArrayRead(patch->localY, &localY);CHKERRQ(ierr);
1693   ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr);
1694   ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr);
1695   ierr = VecRestoreArrayRead(patch->localY, &localY);CHKERRQ(ierr);
1696 
1697   /* Now we need to send the global BC values through */
1698   ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr);
1699   ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr);
1700   ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
1701   ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr);
1702   for (bc = 0; bc < numBcs; ++bc) {
1703     const PetscInt idx = bcNodes[bc];
1704     if (idx < n) globalY[idx] = globalX[idx];
1705   }
1706 
1707   ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
1708   ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr);
1709   ierr = VecRestoreArray(y, &globalY);CHKERRQ(ierr);
1710 
1711   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
1712   ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
1713   PetscFunctionReturn(0);
1714 }
1715 
1716 static PetscErrorCode PCReset_PATCH(PC pc)
1717 {
1718   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1719   PetscInt       i;
1720   PetscErrorCode ierr;
1721 
1722   PetscFunctionBegin;
1723   /* TODO: Get rid of all these ifs */
1724   ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr);
1725   ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr);
1726   ierr = PetscSectionDestroy(&patch->pointCounts);CHKERRQ(ierr);
1727   ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr);
1728   ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr);
1729   ierr = PetscSectionDestroy(&patch->bcCounts);CHKERRQ(ierr);
1730   ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr);
1731   ierr = ISDestroy(&patch->cells);CHKERRQ(ierr);
1732   ierr = ISDestroy(&patch->points);CHKERRQ(ierr);
1733   ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr);
1734   ierr = ISDestroy(&patch->offs);CHKERRQ(ierr);
1735   ierr = PetscSectionDestroy(&patch->patchSection);CHKERRQ(ierr);
1736   ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr);
1737   ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr);
1738 
1739   if (patch->dofSection) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);}
1740   ierr = PetscFree(patch->dofSection);CHKERRQ(ierr);
1741   ierr = PetscFree(patch->bs);CHKERRQ(ierr);
1742   ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr);
1743   if (patch->cellNodeMap) for (i = 0; i < patch->nsubspaces; i++) {ierr = PetscFree(patch->cellNodeMap[i]);CHKERRQ(ierr);}
1744   ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr);
1745   ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr);
1746 
1747   if (patch->bcs) {
1748     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->bcs[i]);CHKERRQ(ierr);}
1749     ierr = PetscFree(patch->bcs);CHKERRQ(ierr);
1750   }
1751   if (patch->multBcs) {
1752     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->multBcs[i]);CHKERRQ(ierr);}
1753     ierr = PetscFree(patch->multBcs);CHKERRQ(ierr);
1754   }
1755   if (patch->ksp) {
1756     for (i = 0; i < patch->npatch; ++i) {ierr = KSPReset(patch->ksp[i]);CHKERRQ(ierr);}
1757   }
1758 
1759   ierr = VecDestroy(&patch->localX);CHKERRQ(ierr);
1760   ierr = VecDestroy(&patch->localY);CHKERRQ(ierr);
1761   if (patch->patchX) {
1762     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchX[i]);CHKERRQ(ierr);}
1763     ierr = PetscFree(patch->patchX);CHKERRQ(ierr);
1764   }
1765   if (patch->patchY) {
1766     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patchY[i]);CHKERRQ(ierr);}
1767     ierr = PetscFree(patch->patchY);CHKERRQ(ierr);
1768   }
1769   ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr);
1770   if (patch->patch_dof_weights) {
1771     for (i = 0; i < patch->npatch; ++i) {ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);}
1772     ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr);
1773   }
1774   if (patch->mat) {
1775     for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);}
1776     ierr = PetscFree(patch->mat);CHKERRQ(ierr);
1777   }
1778   if (patch->multMat) {
1779     for (i = 0; i < patch->npatch; ++i) {ierr = MatDestroy(&patch->multMat[i]);CHKERRQ(ierr);}
1780     ierr = PetscFree(patch->multMat);CHKERRQ(ierr);
1781   }
1782   ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);
1783   if (patch->userIS) {
1784     for (i = 0; i < patch->npatch; ++i) {ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);}
1785     ierr = PetscFree(patch->userIS);CHKERRQ(ierr);
1786   }
1787   patch->bs          = 0;
1788   patch->cellNodeMap = NULL;
1789   ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr);
1790 
1791   ierr = PetscViewerDestroy(&patch->viewerSection);CHKERRQ(ierr);
1792   PetscFunctionReturn(0);
1793 }
1794 
1795 static PetscErrorCode PCDestroy_PATCH(PC pc)
1796 {
1797   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1798   PetscInt       i;
1799   PetscErrorCode ierr;
1800 
1801   PetscFunctionBegin;
1802   ierr = PCReset_PATCH(pc);CHKERRQ(ierr);
1803   if (patch->ksp) {
1804     for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy(&patch->ksp[i]);CHKERRQ(ierr);}
1805     ierr = PetscFree(patch->ksp);CHKERRQ(ierr);
1806   }
1807   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1808   PetscFunctionReturn(0);
1809 }
1810 
1811 static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc)
1812 {
1813   PC_PATCH            *patch = (PC_PATCH *) pc->data;
1814   PCPatchConstructType patchConstructionType = PC_PATCH_STAR;
1815   char                 sub_mat_type[PETSC_MAX_PATH_LEN];
1816   const char          *prefix;
1817   PetscBool            flg, dimflg, codimflg;
1818   MPI_Comm             comm;
1819   PetscErrorCode       ierr;
1820 
1821   PetscFunctionBegin;
1822   ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr);
1823   ierr = PetscObjectGetOptionsPrefix((PetscObject) pc, &prefix);CHKERRQ(ierr);
1824   ierr = PetscOptionsHead(PetscOptionsObject, "Vertex-patch Additive Schwarz options");CHKERRQ(ierr);
1825   ierr = PetscOptionsBool("-pc_patch_save_operators",  "Store all patch operators for lifetime of PC?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr);
1826   ierr = PetscOptionsBool("-pc_patch_partition_of_unity", "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr);
1827   ierr = PetscOptionsBool("-pc_patch_multiplicative", "Gauss-Seidel instead of Jacobi?", "PCPatchSetMultiplicative", patch->multiplicative, &patch->multiplicative, &flg);CHKERRQ(ierr);
1828   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);
1829   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);
1830   if (dimflg && codimflg) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension");CHKERRQ(ierr);
1831   ierr = PetscOptionsEnum("-pc_patch_construct_type", "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr);
1832   if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);}
1833   ierr = PetscOptionsInt("-pc_patch_vanka_dim", "Topological dimension of entities for Vanka to ignore", "PCPATCH", patch->vankadim, &patch->vankadim, &flg);CHKERRQ(ierr);
1834   ierr = PetscOptionsInt("-pc_patch_ignore_dim", "Topological dimension of entities for completion to ignore", "PCPATCH", patch->ignoredim, &patch->ignoredim, &flg);CHKERRQ(ierr);
1835   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);
1836   if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);}
1837   ierr = PetscOptionsBool("-pc_patch_symmetrise_sweep", "Go start->end, end->start?", "PCPATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr);
1838   ierr = PetscOptionsInt("-pc_patch_exclude_subspace", "What subspace (if any) to exclude in construction?", "PCPATCH", patch->exclude_subspace, &patch->exclude_subspace, &flg);CHKERRQ(ierr);
1839 
1840   ierr = PetscOptionsBool("-pc_patch_patches_view", "Print out information during patch construction", "PCPATCH", patch->viewPatches, &patch->viewPatches, &flg);CHKERRQ(ierr);
1841   ierr = PetscOptionsGetViewer(comm, prefix, "-pc_patch_cells_view",   &patch->viewerCells,   &patch->formatCells,   &patch->viewCells);CHKERRQ(ierr);
1842   ierr = PetscOptionsGetViewer(comm, prefix, "-pc_patch_points_view",  &patch->viewerPoints,  &patch->formatPoints,  &patch->viewPoints);CHKERRQ(ierr);
1843   ierr = PetscOptionsGetViewer(comm, prefix, "-pc_patch_section_view", &patch->viewerSection, &patch->formatSection, &patch->viewSection);CHKERRQ(ierr);
1844   ierr = PetscOptionsGetViewer(comm, prefix, "-pc_patch_sub_mat_view", &patch->viewerMatrix,  &patch->formatMatrix,  &patch->viewMatrix);CHKERRQ(ierr);
1845   ierr = PetscOptionsTail();CHKERRQ(ierr);
1846   patch->optionsSet = PETSC_TRUE;
1847   PetscFunctionReturn(0);
1848 }
1849 
1850 static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc)
1851 {
1852   PC_PATCH          *patch = (PC_PATCH*) pc->data;
1853   KSPConvergedReason reason;
1854   PetscInt           i;
1855   PetscErrorCode     ierr;
1856 
1857   PetscFunctionBegin;
1858   for (i = 0; i < patch->npatch; ++i) {
1859     ierr = KSPSetUp(patch->ksp[i]);CHKERRQ(ierr);
1860     ierr = KSPGetConvergedReason(patch->ksp[i], &reason);CHKERRQ(ierr);
1861     if (reason == KSP_DIVERGED_PCSETUP_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1862   }
1863   PetscFunctionReturn(0);
1864 }
1865 
1866 static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer)
1867 {
1868   PC_PATCH      *patch = (PC_PATCH *) pc->data;
1869   PetscViewer    sviewer;
1870   PetscBool      isascii;
1871   PetscMPIInt    rank;
1872   PetscErrorCode ierr;
1873 
1874   PetscFunctionBegin;
1875   /* TODO Redo tabbing with set tbas in new style */
1876   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
1877   if (!isascii) PetscFunctionReturn(0);
1878   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr);
1879   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1880   ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr);
1881   if (patch->multiplicative) {ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: multiplicative\n");CHKERRQ(ierr);}
1882   else                       {ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr);}
1883   if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);}
1884   else                           {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);}
1885   if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);}
1886   else                         {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);}
1887   if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);}
1888   else                        {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);}
1889   if (patch->patchconstructop == PCPatchConstruct_Star)       {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);}
1890   else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);}
1891   else if (patch->patchconstructop == PCPatchConstruct_User)  {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);}
1892   else                                                        {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);}
1893   ierr = PetscViewerASCIIPrintf(viewer, "KSP on patches (all same):\n");CHKERRQ(ierr);
1894   if (patch->ksp) {
1895     ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
1896     if (!rank) {
1897       ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr);
1898       ierr = KSPView(patch->ksp[0], sviewer);CHKERRQ(ierr);
1899       ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr);
1900     }
1901     ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
1902   } else {
1903     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1904     ierr = PetscViewerASCIIPrintf(viewer, "KSP not yet set.\n");CHKERRQ(ierr);
1905     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1906   }
1907   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1908   PetscFunctionReturn(0);
1909 }
1910 
1911 /*MC
1912   PCPATCH = "patch" - A PC object that encapsulates flexible definition of blocks for overlapping and non-overlapping
1913                       small block additive and multiplicative preconditioners. Block definition is based on topology from
1914                       a DM and equation numbering from a PetscSection.
1915 
1916   Options Database Keys:
1917 + -pc_patch_cells_view   - Views the process local cell numbers for each patch
1918 . -pc_patch_points_view  - Views the process local mesh point numbers for each patch
1919 . -pc_patch_g2l_view     - Views the map between global dofs and patch local dofs for each patch
1920 . -pc_patch_patches_view - Views the global dofs associated with each patch and its boundary
1921 - -pc_patch_sub_mat_view - Views the matrix associated with each patch
1922 
1923   Level: intermediate
1924 
1925 .seealso: PCType, PCCreate(), PCSetType()
1926 M*/
1927 PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc)
1928 {
1929   PC_PATCH      *patch;
1930   PetscErrorCode ierr;
1931 
1932   PetscFunctionBegin;
1933   ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr);
1934 
1935   /* Set some defaults */
1936   patch->combined           = PETSC_FALSE;
1937   patch->save_operators     = PETSC_TRUE;
1938   patch->partition_of_unity = PETSC_FALSE;
1939   patch->multiplicative     = PETSC_FALSE;
1940   patch->codim              = -1;
1941   patch->dim                = -1;
1942   patch->exclude_subspace   = -1;
1943   patch->vankadim           = -1;
1944   patch->ignoredim          = -1;
1945   patch->patchconstructop   = PCPatchConstruct_Star;
1946   patch->symmetrise_sweep   = PETSC_FALSE;
1947   patch->npatch             = 0;
1948   patch->userIS             = NULL;
1949   patch->optionsSet         = PETSC_FALSE;
1950   patch->iterationSet       = NULL;
1951   patch->user_patches       = PETSC_FALSE;
1952   ierr = PetscStrallocpy(MATDENSE, (char **) &patch->sub_mat_type);CHKERRQ(ierr);
1953   patch->viewPatches        = PETSC_FALSE;
1954   patch->viewCells          = PETSC_FALSE;
1955   patch->viewPoints         = PETSC_FALSE;
1956   patch->viewSection        = PETSC_FALSE;
1957   patch->viewMatrix         = PETSC_FALSE;
1958 
1959   pc->data                 = (void *) patch;
1960   pc->ops->apply           = PCApply_PATCH;
1961   pc->ops->applytranspose  = 0; /* PCApplyTranspose_PATCH; */
1962   pc->ops->setup           = PCSetUp_PATCH;
1963   pc->ops->reset           = PCReset_PATCH;
1964   pc->ops->destroy         = PCDestroy_PATCH;
1965   pc->ops->setfromoptions  = PCSetFromOptions_PATCH;
1966   pc->ops->setuponblocks   = PCSetUpOnBlocks_PATCH;
1967   pc->ops->view            = PCView_PATCH;
1968   pc->ops->applyrichardson = 0;
1969 
1970   PetscFunctionReturn(0);
1971 }
1972