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