xref: /petsc/src/ksp/pc/impls/patch/pcpatch.c (revision 642283e9ad0dfe8629331f3ed3cf8d5bb073eab4)
14bbf5ea8SMatthew G. Knepley #include <petsc/private/pcpatchimpl.h>     /*I "petscpc.h" I*/
24bbf5ea8SMatthew G. Knepley #include <petscdmplex.h>
34bbf5ea8SMatthew G. Knepley #include <petscsf.h>
44bbf5ea8SMatthew G. Knepley #include <petscbt.h>
54bbf5ea8SMatthew G. Knepley 
64bbf5ea8SMatthew G. Knepley PetscLogEvent PC_Patch_CreatePatches, PC_Patch_ComputeOp, PC_Patch_Solve, PC_Patch_Scatter, PC_Patch_Apply, PC_Patch_Prealloc;
74bbf5ea8SMatthew G. Knepley 
84bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchConstruct_Star(void *vpatch, DM dm, PetscInt point, PetscHashI ht)
94bbf5ea8SMatthew G. Knepley {
104bbf5ea8SMatthew G. Knepley   PetscInt       starSize;
114bbf5ea8SMatthew G. Knepley   PetscInt      *star = NULL, si;
124bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
134bbf5ea8SMatthew G. Knepley 
144bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
154bbf5ea8SMatthew G. Knepley   PetscHashIClear(ht);
164bbf5ea8SMatthew G. Knepley   /* To start with, add the point we care about */
174bbf5ea8SMatthew G. Knepley   PetscHashIAdd(ht, point, 0);
184bbf5ea8SMatthew G. Knepley   /* Loop over all the points that this point connects to */
194bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
204bbf5ea8SMatthew G. Knepley   for (si = 0; si < starSize; si += 2) {PetscHashIAdd(ht, star[si], 0);}
214bbf5ea8SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
224bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
234bbf5ea8SMatthew G. Knepley }
244bbf5ea8SMatthew G. Knepley 
254bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchConstruct_Vanka(void *vpatch, DM dm, PetscInt point, PetscHashI ht)
264bbf5ea8SMatthew G. Knepley {
274bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) vpatch;
284bbf5ea8SMatthew G. Knepley   PetscInt       starSize;
294bbf5ea8SMatthew G. Knepley   PetscInt      *star = NULL;
304bbf5ea8SMatthew G. Knepley   PetscBool      shouldIgnore = PETSC_FALSE;
314bbf5ea8SMatthew G. Knepley   PetscInt       cStart, cEnd, iStart, iEnd, si;
324bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
334bbf5ea8SMatthew G. Knepley 
344bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
354bbf5ea8SMatthew G. Knepley   PetscHashIClear(ht);
364bbf5ea8SMatthew G. Knepley   /* To start with, add the point we care about */
374bbf5ea8SMatthew G. Knepley   PetscHashIAdd(ht, point, 0);
384bbf5ea8SMatthew G. Knepley   /* Should we ignore any points of a certain dimension? */
394bbf5ea8SMatthew G. Knepley   if (patch->vankadim >= 0) {
404bbf5ea8SMatthew G. Knepley     shouldIgnore = PETSC_TRUE;
414bbf5ea8SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, patch->vankadim, &iStart, &iEnd);CHKERRQ(ierr);
424bbf5ea8SMatthew G. Knepley   }
434bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
444bbf5ea8SMatthew G. Knepley   /* Loop over all the cells that this point connects to */
454bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
464bbf5ea8SMatthew G. Knepley   for (si = 0; si < starSize; si += 2) {
474bbf5ea8SMatthew G. Knepley     const PetscInt cell = star[si];
484bbf5ea8SMatthew G. Knepley     PetscInt       closureSize;
494bbf5ea8SMatthew G. Knepley     PetscInt      *closure = NULL, ci;
504bbf5ea8SMatthew G. Knepley 
514bbf5ea8SMatthew G. Knepley     if (cell < cStart || cell >= cEnd) continue;
524bbf5ea8SMatthew G. Knepley     /* now loop over all entities in the closure of that cell */
534bbf5ea8SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
544bbf5ea8SMatthew G. Knepley     for (ci = 0; ci < closureSize; ci += 2) {
554bbf5ea8SMatthew G. Knepley       const PetscInt newpoint = closure[ci];
564bbf5ea8SMatthew G. Knepley 
574bbf5ea8SMatthew G. Knepley       /* We've been told to ignore entities of this type.*/
584bbf5ea8SMatthew G. Knepley       if (shouldIgnore && newpoint >= iStart && newpoint < iEnd) continue;
594bbf5ea8SMatthew G. Knepley       PetscHashIAdd(ht, newpoint, 0);
604bbf5ea8SMatthew G. Knepley     }
614bbf5ea8SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, cell, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
624bbf5ea8SMatthew G. Knepley   }
634bbf5ea8SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
644bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
654bbf5ea8SMatthew G. Knepley }
664bbf5ea8SMatthew G. Knepley 
674bbf5ea8SMatthew G. Knepley /* The user's already set the patches in patch->userIS. Build the hash tables */
684bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchConstruct_User(void *vpatch, DM dm, PetscInt point, PetscHashI ht)
694bbf5ea8SMatthew G. Knepley {
704bbf5ea8SMatthew G. Knepley   PC_PATCH       *patch   = (PC_PATCH *) vpatch;
714bbf5ea8SMatthew G. Knepley   IS              patchis = patch->userIS[point];
724bbf5ea8SMatthew G. Knepley   PetscInt        n;
734bbf5ea8SMatthew G. Knepley   const PetscInt *patchdata;
744bbf5ea8SMatthew G. Knepley   PetscInt        pStart, pEnd, i;
754bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
764bbf5ea8SMatthew G. Knepley 
774bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
784bbf5ea8SMatthew G. Knepley   PetscHashIClear(ht);
794bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);
804bbf5ea8SMatthew G. Knepley   ierr = ISGetLocalSize(patchis, &n);CHKERRQ(ierr);
814bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(patchis, &patchdata);CHKERRQ(ierr);
824bbf5ea8SMatthew G. Knepley   for (i = 0; i < n; ++i) {
834bbf5ea8SMatthew G. Knepley     const PetscInt ownedpoint = patchdata[i];
844bbf5ea8SMatthew G. Knepley 
854bbf5ea8SMatthew G. Knepley     if (ownedpoint < pStart || ownedpoint >= pEnd) {
864bbf5ea8SMatthew G. Knepley       SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Mesh point %D was not in [%D, %D)", ownedpoint, pStart, pEnd);
874bbf5ea8SMatthew G. Knepley     }
884bbf5ea8SMatthew G. Knepley     PetscHashIAdd(ht, ownedpoint, 0);
894bbf5ea8SMatthew G. Knepley   }
904bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(patchis, &patchdata);CHKERRQ(ierr);
914bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
924bbf5ea8SMatthew G. Knepley }
934bbf5ea8SMatthew G. Knepley 
944bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchCreateDefaultSF_Private(PC pc, PetscInt n, const PetscSF *sf, const PetscInt *bs)
954bbf5ea8SMatthew G. Knepley {
964bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
974bbf5ea8SMatthew G. Knepley   PetscInt       i;
984bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
994bbf5ea8SMatthew G. Knepley 
1004bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
1014bbf5ea8SMatthew G. Knepley   if (n == 1 && bs[0] == 1) {
1024bbf5ea8SMatthew G. Knepley     patch->defaultSF = sf[0];
1034bbf5ea8SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) patch->defaultSF);CHKERRQ(ierr);
1044bbf5ea8SMatthew G. Knepley   } else {
1054bbf5ea8SMatthew G. Knepley     PetscInt     allRoots = 0, allLeaves = 0;
1064bbf5ea8SMatthew G. Knepley     PetscInt     leafOffset = 0;
1074bbf5ea8SMatthew G. Knepley     PetscInt    *ilocal = NULL;
1084bbf5ea8SMatthew G. Knepley     PetscSFNode *iremote = NULL;
1094bbf5ea8SMatthew G. Knepley     PetscInt    *remoteOffsets = NULL;
1104bbf5ea8SMatthew G. Knepley     PetscInt     index = 0;
1114bbf5ea8SMatthew G. Knepley     PetscHashI   rankToIndex;
1124bbf5ea8SMatthew G. Knepley     PetscInt     numRanks = 0;
1134bbf5ea8SMatthew G. Knepley     PetscSFNode *remote = NULL;
1144bbf5ea8SMatthew G. Knepley     PetscSF      rankSF;
1154bbf5ea8SMatthew G. Knepley     PetscInt    *ranks = NULL;
1164bbf5ea8SMatthew G. Knepley     PetscInt    *offsets = NULL;
1174bbf5ea8SMatthew G. Knepley     MPI_Datatype contig;
1184bbf5ea8SMatthew G. Knepley     PetscHashI   ranksUniq;
1194bbf5ea8SMatthew G. Knepley 
1204bbf5ea8SMatthew G. Knepley     /* First figure out how many dofs there are in the concatenated numbering.
1214bbf5ea8SMatthew G. Knepley      * allRoots: number of owned global dofs;
1224bbf5ea8SMatthew G. Knepley      * allLeaves: number of visible dofs (global + ghosted).
1234bbf5ea8SMatthew G. Knepley      */
1244bbf5ea8SMatthew G. Knepley     for (i = 0; i < n; ++i) {
1254bbf5ea8SMatthew G. Knepley       PetscInt nroots, nleaves;
1264bbf5ea8SMatthew G. Knepley 
1274bbf5ea8SMatthew G. Knepley       ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, NULL, NULL);CHKERRQ(ierr);
1284bbf5ea8SMatthew G. Knepley       allRoots  += nroots * bs[i];
1294bbf5ea8SMatthew G. Knepley       allLeaves += nleaves * bs[i];
1304bbf5ea8SMatthew G. Knepley     }
1314bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(allLeaves, &ilocal);CHKERRQ(ierr);
1324bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(allLeaves, &iremote);CHKERRQ(ierr);
1334bbf5ea8SMatthew G. Knepley     /* Now build an SF that just contains process connectivity. */
1344bbf5ea8SMatthew G. Knepley     PetscHashICreate(ranksUniq);
1354bbf5ea8SMatthew G. Knepley     for (i = 0; i < n; ++i) {
1364bbf5ea8SMatthew G. Knepley       const PetscMPIInt *ranks = NULL;
1374bbf5ea8SMatthew G. Knepley       PetscInt           nranks, j;
1384bbf5ea8SMatthew G. Knepley 
1394bbf5ea8SMatthew G. Knepley       ierr = PetscSFSetUp(sf[i]);CHKERRQ(ierr);
1404bbf5ea8SMatthew G. Knepley       ierr = PetscSFGetRanks(sf[i], &nranks, &ranks, NULL, NULL, NULL);CHKERRQ(ierr);
1414bbf5ea8SMatthew G. Knepley       /* These are all the ranks who communicate with me. */
1424bbf5ea8SMatthew G. Knepley       for (j = 0; j < nranks; ++j) {
1434bbf5ea8SMatthew G. Knepley         PetscHashIAdd(ranksUniq, (PetscInt) ranks[j], 0);
1444bbf5ea8SMatthew G. Knepley       }
1454bbf5ea8SMatthew G. Knepley     }
1464bbf5ea8SMatthew G. Knepley     PetscHashISize(ranksUniq, numRanks);
1474bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(numRanks, &remote);CHKERRQ(ierr);
1484bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(numRanks, &ranks);CHKERRQ(ierr);
1494bbf5ea8SMatthew G. Knepley     PetscHashIGetKeys(ranksUniq, &index, ranks);
1504bbf5ea8SMatthew G. Knepley 
1514bbf5ea8SMatthew G. Knepley     PetscHashICreate(rankToIndex);
1524bbf5ea8SMatthew G. Knepley     for (i = 0; i < numRanks; ++i) {
1534bbf5ea8SMatthew G. Knepley       remote[i].rank  = ranks[i];
1544bbf5ea8SMatthew G. Knepley       remote[i].index = 0;
1554bbf5ea8SMatthew G. Knepley       PetscHashIAdd(rankToIndex, ranks[i], i);
1564bbf5ea8SMatthew G. Knepley     }
1574bbf5ea8SMatthew G. Knepley     ierr = PetscFree(ranks);CHKERRQ(ierr);
1584bbf5ea8SMatthew G. Knepley     PetscHashIDestroy(ranksUniq);
1594bbf5ea8SMatthew G. Knepley     ierr = PetscSFCreate(PetscObjectComm((PetscObject) pc), &rankSF);CHKERRQ(ierr);
1604bbf5ea8SMatthew G. Knepley     ierr = PetscSFSetGraph(rankSF, 1, numRanks, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1614bbf5ea8SMatthew G. Knepley     ierr = PetscSFSetUp(rankSF);CHKERRQ(ierr);
1624bbf5ea8SMatthew G. Knepley     /* OK, use it to communicate the root offset on the remote
1634bbf5ea8SMatthew G. Knepley      * processes for each subspace. */
1644bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(n, &offsets);CHKERRQ(ierr);
1654bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(n*numRanks, &remoteOffsets);CHKERRQ(ierr);
1664bbf5ea8SMatthew G. Knepley 
1674bbf5ea8SMatthew G. Knepley     offsets[0] = 0;
1684bbf5ea8SMatthew G. Knepley     for (i = 1; i < n; ++i) {
1694bbf5ea8SMatthew G. Knepley       PetscInt nroots;
1704bbf5ea8SMatthew G. Knepley 
1714bbf5ea8SMatthew G. Knepley       ierr = PetscSFGetGraph(sf[i-1], &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1724bbf5ea8SMatthew G. Knepley       offsets[i] = offsets[i-1] + nroots*bs[i-1];
1734bbf5ea8SMatthew G. Knepley     }
1744bbf5ea8SMatthew G. Knepley     /* Offsets are the offsets on the current process of the
1754bbf5ea8SMatthew G. Knepley      * global dof numbering for the subspaces. */
1764bbf5ea8SMatthew G. Knepley     ierr = MPI_Type_contiguous(n, MPIU_INT, &contig);CHKERRQ(ierr);
1774bbf5ea8SMatthew G. Knepley     ierr = MPI_Type_commit(&contig);CHKERRQ(ierr);
1784bbf5ea8SMatthew G. Knepley 
1794bbf5ea8SMatthew G. Knepley     ierr = PetscSFBcastBegin(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr);
1804bbf5ea8SMatthew G. Knepley     ierr = PetscSFBcastEnd(rankSF, contig, offsets, remoteOffsets);CHKERRQ(ierr);
1814bbf5ea8SMatthew G. Knepley     ierr = MPI_Type_free(&contig);CHKERRQ(ierr);
1824bbf5ea8SMatthew G. Knepley     ierr = PetscFree(offsets);CHKERRQ(ierr);
1834bbf5ea8SMatthew G. Knepley     ierr = PetscSFDestroy(&rankSF);CHKERRQ(ierr);
1844bbf5ea8SMatthew G. Knepley     /* Now remoteOffsets contains the offsets on the remote
1854bbf5ea8SMatthew G. Knepley      * processes who communicate with me.  So now we can
1864bbf5ea8SMatthew G. Knepley      * concatenate the list of SFs into a single one. */
1874bbf5ea8SMatthew G. Knepley     index = 0;
1884bbf5ea8SMatthew G. Knepley     for (i = 0; i < n; ++i) {
1894bbf5ea8SMatthew G. Knepley       const PetscSFNode *remote = NULL;
1904bbf5ea8SMatthew G. Knepley       const PetscInt    *local  = NULL;
1914bbf5ea8SMatthew G. Knepley       PetscInt           nroots, nleaves, j;
1924bbf5ea8SMatthew G. Knepley 
1934bbf5ea8SMatthew G. Knepley       ierr = PetscSFGetGraph(sf[i], &nroots, &nleaves, &local, &remote);CHKERRQ(ierr);
1944bbf5ea8SMatthew G. Knepley       for (j = 0; j < nleaves; ++j) {
1954bbf5ea8SMatthew G. Knepley         PetscInt rank = remote[j].rank;
1964bbf5ea8SMatthew G. Knepley         PetscInt idx, rootOffset, k;
1974bbf5ea8SMatthew G. Knepley 
1984bbf5ea8SMatthew G. Knepley         PetscHashIMap(rankToIndex, rank, idx);
1994bbf5ea8SMatthew G. Knepley         if (idx == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Didn't find rank, huh?");
2004bbf5ea8SMatthew G. Knepley         /* Offset on given rank for ith subspace */
2014bbf5ea8SMatthew G. Knepley         rootOffset = remoteOffsets[n*idx + i];
2024bbf5ea8SMatthew G. Knepley         for (k = 0; k < bs[i]; ++k) {
2034bbf5ea8SMatthew G. Knepley           ilocal[index]        = local[j]*bs[i] + k + leafOffset;
2044bbf5ea8SMatthew G. Knepley           iremote[index].rank  = remote[j].rank;
2054bbf5ea8SMatthew G. Knepley           iremote[index].index = remote[j].index*bs[i] + k + rootOffset;
2064bbf5ea8SMatthew G. Knepley           ++index;
2074bbf5ea8SMatthew G. Knepley         }
2084bbf5ea8SMatthew G. Knepley       }
2094bbf5ea8SMatthew G. Knepley       leafOffset += nleaves * bs[i];
2104bbf5ea8SMatthew G. Knepley     }
2114bbf5ea8SMatthew G. Knepley     PetscHashIDestroy(rankToIndex);
2124bbf5ea8SMatthew G. Knepley     ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
2134bbf5ea8SMatthew G. Knepley     ierr = PetscSFCreate(PetscObjectComm((PetscObject)pc), &patch->defaultSF);CHKERRQ(ierr);
2144bbf5ea8SMatthew G. Knepley     ierr = PetscSFSetGraph(patch->defaultSF, allRoots, allLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
2154bbf5ea8SMatthew G. Knepley   }
2164bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2174bbf5ea8SMatthew G. Knepley }
2184bbf5ea8SMatthew G. Knepley 
2194bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2204bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetSaveOperators(PC pc, PetscBool flg)
2214bbf5ea8SMatthew G. Knepley {
2224bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2234bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2244bbf5ea8SMatthew G. Knepley   patch->save_operators = flg;
2254bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2264bbf5ea8SMatthew G. Knepley }
2274bbf5ea8SMatthew G. Knepley 
2284bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2294bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetSaveOperators(PC pc, PetscBool *flg)
2304bbf5ea8SMatthew G. Knepley {
2314bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2324bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2334bbf5ea8SMatthew G. Knepley   *flg = patch->save_operators;
2344bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2354bbf5ea8SMatthew G. Knepley }
2364bbf5ea8SMatthew G. Knepley 
2374bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2384bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetPartitionOfUnity(PC pc, PetscBool flg)
2394bbf5ea8SMatthew G. Knepley {
2404bbf5ea8SMatthew G. Knepley     PC_PATCH *patch = (PC_PATCH *) pc->data;
2414bbf5ea8SMatthew G. Knepley     PetscFunctionBegin;
2424bbf5ea8SMatthew G. Knepley     patch->partition_of_unity = flg;
2434bbf5ea8SMatthew G. Knepley     PetscFunctionReturn(0);
2444bbf5ea8SMatthew G. Knepley }
2454bbf5ea8SMatthew G. Knepley 
2464bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2474bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetPartitionOfUnity(PC pc, PetscBool *flg)
2484bbf5ea8SMatthew G. Knepley {
2494bbf5ea8SMatthew G. Knepley     PC_PATCH *patch = (PC_PATCH *) pc->data;
2504bbf5ea8SMatthew G. Knepley     PetscFunctionBegin;
2514bbf5ea8SMatthew G. Knepley     *flg = patch->partition_of_unity;
2524bbf5ea8SMatthew G. Knepley     PetscFunctionReturn(0);
2534bbf5ea8SMatthew G. Knepley }
2544bbf5ea8SMatthew G. Knepley 
2554bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2564bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetMultiplicative(PC pc, PetscBool flg)
2574bbf5ea8SMatthew G. Knepley {
2584bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2594bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2604bbf5ea8SMatthew G. Knepley   patch->multiplicative = flg;
2614bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2624bbf5ea8SMatthew G. Knepley }
2634bbf5ea8SMatthew G. Knepley 
2644bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2654bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetMultiplicative(PC pc, PetscBool *flg)
2664bbf5ea8SMatthew G. Knepley {
2674bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2684bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2694bbf5ea8SMatthew G. Knepley   *flg = patch->multiplicative;
2704bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2714bbf5ea8SMatthew G. Knepley }
2724bbf5ea8SMatthew G. Knepley 
2734bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2744bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetSubMatType(PC pc, MatType sub_mat_type)
2754bbf5ea8SMatthew G. Knepley {
2764bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2774bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
2784bbf5ea8SMatthew G. Knepley 
2794bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2804bbf5ea8SMatthew G. Knepley   if (patch->sub_mat_type) {ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);}
2814bbf5ea8SMatthew G. Knepley   ierr = PetscStrallocpy(sub_mat_type, (char **) &patch->sub_mat_type);CHKERRQ(ierr);
2824bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2834bbf5ea8SMatthew G. Knepley }
2844bbf5ea8SMatthew G. Knepley 
2854bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2864bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetSubMatType(PC pc, MatType *sub_mat_type)
2874bbf5ea8SMatthew G. Knepley {
2884bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
2894bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
2904bbf5ea8SMatthew G. Knepley   *sub_mat_type = patch->sub_mat_type;
2914bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
2924bbf5ea8SMatthew G. Knepley }
2934bbf5ea8SMatthew G. Knepley 
2944bbf5ea8SMatthew G. Knepley /* TODO: Docs */
2954bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetCellNumbering(PC pc, PetscSection cellNumbering)
2964bbf5ea8SMatthew G. Knepley {
2974bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
2984bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
2994bbf5ea8SMatthew G. Knepley 
3004bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3014bbf5ea8SMatthew G. Knepley   patch->cellNumbering = cellNumbering;
3024bbf5ea8SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) cellNumbering);CHKERRQ(ierr);
3034bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3044bbf5ea8SMatthew G. Knepley }
3054bbf5ea8SMatthew G. Knepley 
3064bbf5ea8SMatthew G. Knepley /* TODO: Docs */
3074bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetCellNumbering(PC pc, PetscSection *cellNumbering)
3084bbf5ea8SMatthew G. Knepley {
3094bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
3104bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3114bbf5ea8SMatthew G. Knepley   *cellNumbering = patch->cellNumbering;
3124bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3134bbf5ea8SMatthew G. Knepley }
3144bbf5ea8SMatthew G. Knepley 
3154bbf5ea8SMatthew G. Knepley /* TODO: Docs */
3164bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetConstructType(PC pc, PCPatchConstructType ctype, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *ctx)
3174bbf5ea8SMatthew G. Knepley {
3184bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
3194bbf5ea8SMatthew G. Knepley 
3204bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3214bbf5ea8SMatthew G. Knepley   patch->ctype = ctype;
3224bbf5ea8SMatthew G. Knepley   switch (ctype) {
3234bbf5ea8SMatthew G. Knepley   case PC_PATCH_STAR:
3244bbf5ea8SMatthew G. Knepley     patch->patchconstructop = PCPatchConstruct_Star;
3254bbf5ea8SMatthew G. Knepley     break;
3264bbf5ea8SMatthew G. Knepley   case PC_PATCH_VANKA:
3274bbf5ea8SMatthew G. Knepley     patch->patchconstructop = PCPatchConstruct_Vanka;
3284bbf5ea8SMatthew G. Knepley     break;
3294bbf5ea8SMatthew G. Knepley   case PC_PATCH_USER:
3304bbf5ea8SMatthew G. Knepley   case PC_PATCH_PYTHON:
3314bbf5ea8SMatthew G. Knepley     patch->user_patches            = PETSC_TRUE;
3324bbf5ea8SMatthew G. Knepley     patch->patchconstructop        = PCPatchConstruct_User;
3334bbf5ea8SMatthew G. Knepley     patch->userpatchconstructionop = func;
3344bbf5ea8SMatthew G. Knepley     patch->userpatchconstructctx   = ctx;
3354bbf5ea8SMatthew G. Knepley     break;
3364bbf5ea8SMatthew G. Knepley   default:
3374bbf5ea8SMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype);
3384bbf5ea8SMatthew G. Knepley   }
3394bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3404bbf5ea8SMatthew G. Knepley }
3414bbf5ea8SMatthew G. Knepley 
3424bbf5ea8SMatthew G. Knepley /* TODO: Docs */
3434bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchGetConstructType(PC pc, PCPatchConstructType *ctype, PetscErrorCode (**func)(PC, PetscInt *, IS **, IS *, void *), void **ctx)
3444bbf5ea8SMatthew G. Knepley {
3454bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
3464bbf5ea8SMatthew G. Knepley 
3474bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3484bbf5ea8SMatthew G. Knepley   *ctype = patch->ctype;
3494bbf5ea8SMatthew G. Knepley   switch (patch->ctype) {
3504bbf5ea8SMatthew G. Knepley   case PC_PATCH_STAR:
3514bbf5ea8SMatthew G. Knepley   case PC_PATCH_VANKA:
3524bbf5ea8SMatthew G. Knepley     break;
3534bbf5ea8SMatthew G. Knepley   case PC_PATCH_USER:
3544bbf5ea8SMatthew G. Knepley   case PC_PATCH_PYTHON:
3554bbf5ea8SMatthew G. Knepley     *func = patch->userpatchconstructionop;
3564bbf5ea8SMatthew G. Knepley     *ctx  = patch->userpatchconstructctx;
3574bbf5ea8SMatthew G. Knepley     break;
3584bbf5ea8SMatthew G. Knepley   default:
3594bbf5ea8SMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_USER, "Unknown patch construction type %D", (PetscInt) patch->ctype);
3604bbf5ea8SMatthew G. Knepley   }
3614bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
3624bbf5ea8SMatthew G. Knepley }
3634bbf5ea8SMatthew G. Knepley 
3644bbf5ea8SMatthew G. Knepley /* TODO: Docs */
3654bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetDiscretisationInfo(PC pc, PetscInt nsubspaces, DM *dms, PetscInt *bs, PetscInt *nodesPerCell, const PetscInt **cellNodeMap,
3664bbf5ea8SMatthew G. Knepley                                             const PetscInt *subspaceOffsets, PetscInt numGhostBcs, const PetscInt *ghostBcNodes, PetscInt numGlobalBcs, const PetscInt *globalBcNodes)
3674bbf5ea8SMatthew G. Knepley {
3684bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
3694bbf5ea8SMatthew G. Knepley   PetscSF       *sfs;
3704bbf5ea8SMatthew G. Knepley   PetscInt       i;
3714bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
3724bbf5ea8SMatthew G. Knepley 
3734bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
3744bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &sfs);CHKERRQ(ierr);
3754bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &patch->dofSection);CHKERRQ(ierr);
3764bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &patch->bs);CHKERRQ(ierr);
3774bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &patch->nodesPerCell);CHKERRQ(ierr);
3784bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces, &patch->cellNodeMap);CHKERRQ(ierr);
3794bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(nsubspaces+1, &patch->subspaceOffsets);CHKERRQ(ierr);
3804bbf5ea8SMatthew G. Knepley 
3814bbf5ea8SMatthew G. Knepley   patch->nsubspaces       = nsubspaces;
3824bbf5ea8SMatthew G. Knepley   patch->totalDofsPerCell = 0;
3834bbf5ea8SMatthew G. Knepley   for (i = 0; i < nsubspaces; ++i) {
3844bbf5ea8SMatthew G. Knepley     ierr = DMGetDefaultSection(dms[i], &patch->dofSection[i]);CHKERRQ(ierr);
3854bbf5ea8SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) patch->dofSection[i]);CHKERRQ(ierr);
3864bbf5ea8SMatthew G. Knepley     ierr = DMGetDefaultSF(dms[i], &sfs[i]);CHKERRQ(ierr);
3874bbf5ea8SMatthew G. Knepley     patch->bs[i]              = bs[i];
3884bbf5ea8SMatthew G. Knepley     patch->nodesPerCell[i]    = nodesPerCell[i];
3894bbf5ea8SMatthew G. Knepley     patch->totalDofsPerCell  += nodesPerCell[i]*bs[i];
3904bbf5ea8SMatthew G. Knepley     patch->cellNodeMap[i]     = cellNodeMap[i];
3914bbf5ea8SMatthew G. Knepley     patch->subspaceOffsets[i] = subspaceOffsets[i];
3924bbf5ea8SMatthew G. Knepley   }
3934bbf5ea8SMatthew G. Knepley   ierr = PCPatchCreateDefaultSF_Private(pc, nsubspaces, sfs, patch->bs);CHKERRQ(ierr);
3944bbf5ea8SMatthew G. Knepley   ierr = PetscFree(sfs);CHKERRQ(ierr);
3954bbf5ea8SMatthew G. Knepley 
3964bbf5ea8SMatthew G. Knepley   patch->subspaceOffsets[nsubspaces] = subspaceOffsets[nsubspaces];
3974bbf5ea8SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGhostBcs, ghostBcNodes, PETSC_COPY_VALUES, &patch->ghostBcNodes);CHKERRQ(ierr);
3984bbf5ea8SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalBcs, globalBcNodes, PETSC_COPY_VALUES, &patch->globalBcNodes);CHKERRQ(ierr);
3994bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
4004bbf5ea8SMatthew G. Knepley }
4014bbf5ea8SMatthew G. Knepley 
4024bbf5ea8SMatthew G. Knepley /* TODO: Docs */
4034bbf5ea8SMatthew G. Knepley PetscErrorCode PCPatchSetComputeOperator(PC pc, PetscErrorCode (*func)(PC, Mat, PetscInt, const PetscInt *, PetscInt, const PetscInt *, void *), void *ctx)
4044bbf5ea8SMatthew G. Knepley {
4054bbf5ea8SMatthew G. Knepley   PC_PATCH *patch = (PC_PATCH *) pc->data;
4064bbf5ea8SMatthew G. Knepley 
4074bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
4084bbf5ea8SMatthew G. Knepley   /* User op can assume matrix is zeroed */
4094bbf5ea8SMatthew G. Knepley   patch->usercomputeop  = func;
4104bbf5ea8SMatthew G. Knepley   patch->usercomputectx = ctx;
4114bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
4124bbf5ea8SMatthew G. Knepley }
4134bbf5ea8SMatthew G. Knepley 
4144bbf5ea8SMatthew G. Knepley /* On entry, ht contains the topological entities whose dofs we are responsible for solving for;
4154bbf5ea8SMatthew G. Knepley    on exit, cht contains all the topological entities we need to compute their residuals.
4164bbf5ea8SMatthew G. Knepley    In full generality this should incorporate knowledge of the sparsity pattern of the matrix;
4174bbf5ea8SMatthew G. Knepley    here we assume a standard FE sparsity pattern.*/
4184bbf5ea8SMatthew G. Knepley /* TODO: Use DMPlexGetAdjacency() */
4194bbf5ea8SMatthew G. Knepley /* TODO: Look at temp buffer management for GetClosure() */
4204bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchCompleteCellPatch(DM dm, PetscHashI ht, PetscHashI cht)
4214bbf5ea8SMatthew G. Knepley {
4224bbf5ea8SMatthew G. Knepley   PetscHashIIter hi;
4234bbf5ea8SMatthew G. Knepley   PetscInt       point;
4244bbf5ea8SMatthew G. Knepley   PetscInt      *star = NULL, *closure = NULL;
4254bbf5ea8SMatthew G. Knepley   PetscInt       starSize, closureSize, si, ci;
4264bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
4274bbf5ea8SMatthew G. Knepley 
4284bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
4294bbf5ea8SMatthew G. Knepley   PetscHashIClear(cht);
4304bbf5ea8SMatthew G. Knepley   PetscHashIIterBegin(ht, hi);
4314bbf5ea8SMatthew G. Knepley   while (!PetscHashIIterAtEnd(ht, hi)) {
4324bbf5ea8SMatthew G. Knepley     PetscHashIIterGetKey(ht, hi, point);
4334bbf5ea8SMatthew G. Knepley     PetscHashIIterNext(ht, hi);
4344bbf5ea8SMatthew G. Knepley 
4354bbf5ea8SMatthew G. Knepley     /* Loop over all the cells that this point connects to */
4364bbf5ea8SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_FALSE, &starSize, &star);CHKERRQ(ierr);
4374bbf5ea8SMatthew G. Knepley     for (si = 0; si < starSize; si += 2) {
4384bbf5ea8SMatthew G. Knepley       PetscInt ownedpoint = star[si];
4394bbf5ea8SMatthew G. Knepley       /* now loop over all entities in the closure of that cell */
4404bbf5ea8SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, ownedpoint, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
4414bbf5ea8SMatthew G. Knepley       for (ci = 0; ci < closureSize; ci += 2) {
4424bbf5ea8SMatthew G. Knepley         PetscInt seenpoint = closure[ci];
4434bbf5ea8SMatthew G. Knepley         PetscHashIAdd(cht, seenpoint, 0);
4444bbf5ea8SMatthew G. Knepley       }
4454bbf5ea8SMatthew G. Knepley     }
4464bbf5ea8SMatthew G. Knepley   }
4474bbf5ea8SMatthew G. Knepley   /* Only restore work arrays at very end. */
4484bbf5ea8SMatthew G. Knepley   if (closure) {ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, NULL, &closure);CHKERRQ(ierr);}
4494bbf5ea8SMatthew G. Knepley   if (star)    {ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_FALSE, NULL, &star);CHKERRQ(ierr);}
4504bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
4514bbf5ea8SMatthew G. Knepley }
4524bbf5ea8SMatthew G. Knepley 
4534bbf5ea8SMatthew G. Knepley /* Given a hash table with a set of topological entities (pts), compute the degrees of
4544bbf5ea8SMatthew G. Knepley    freedom in global concatenated numbering on those entities.
4554bbf5ea8SMatthew G. Knepley    For Vanka smoothing, this needs to do something special: ignore dofs of the
4564bbf5ea8SMatthew G. Knepley    constraint subspace on entities that aren't the base entity we're building the patch
4574bbf5ea8SMatthew G. Knepley    around. */
4584bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchGetPointDofs(PC_PATCH *patch, PetscHashI pts, PetscHashI dofs, PetscInt base, PetscInt exclude_subspace)
4594bbf5ea8SMatthew G. Knepley {
4604bbf5ea8SMatthew G. Knepley   PetscHashIIter hi;
4614bbf5ea8SMatthew G. Knepley   PetscInt       ldof, loff;
4624bbf5ea8SMatthew G. Knepley   PetscInt       k, p;
4634bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
4644bbf5ea8SMatthew G. Knepley 
4654bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
4664bbf5ea8SMatthew G. Knepley   PetscHashIClear(dofs);
4674bbf5ea8SMatthew G. Knepley   for (k = 0; k < patch->nsubspaces; ++k) {
4684bbf5ea8SMatthew G. Knepley     PetscSection dofSection     = patch->dofSection[k];
4694bbf5ea8SMatthew G. Knepley     PetscInt     subspaceOffset = patch->subspaceOffsets[k];
4704bbf5ea8SMatthew G. Knepley     PetscInt     bs             = patch->bs[k];
4714bbf5ea8SMatthew G. Knepley     PetscInt     j, l;
4724bbf5ea8SMatthew G. Knepley 
4734bbf5ea8SMatthew G. Knepley     if (k == exclude_subspace) {
4744bbf5ea8SMatthew G. Knepley       /* only get this subspace dofs at the base entity, not any others */
4754bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetDof(dofSection, base, &ldof);CHKERRQ(ierr);
4764bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetOffset(dofSection, base, &loff);CHKERRQ(ierr);
4774bbf5ea8SMatthew G. Knepley       if (0 == ldof) continue;
4784bbf5ea8SMatthew G. Knepley       for (j = loff; j < ldof + loff; ++j) {
4794bbf5ea8SMatthew G. Knepley         for (l = 0; l < bs; ++l) {
4804bbf5ea8SMatthew G. Knepley           PetscInt dof = bs*j + l + subspaceOffset;
4814bbf5ea8SMatthew G. Knepley           PetscHashIAdd(dofs, dof, 0);
4824bbf5ea8SMatthew G. Knepley         }
4834bbf5ea8SMatthew G. Knepley       }
4844bbf5ea8SMatthew G. Knepley       continue; /* skip the other dofs of this subspace */
4854bbf5ea8SMatthew G. Knepley     }
4864bbf5ea8SMatthew G. Knepley 
4874bbf5ea8SMatthew G. Knepley     PetscHashIIterBegin(pts, hi);
4884bbf5ea8SMatthew G. Knepley     while (!PetscHashIIterAtEnd(pts, hi)) {
4894bbf5ea8SMatthew G. Knepley       PetscHashIIterGetKey(pts, hi, p);
4904bbf5ea8SMatthew G. Knepley       PetscHashIIterNext(pts, hi);
4914bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetDof(dofSection, p, &ldof);CHKERRQ(ierr);
4924bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetOffset(dofSection, p, &loff);CHKERRQ(ierr);
4934bbf5ea8SMatthew G. Knepley       if (0 == ldof) continue;
4944bbf5ea8SMatthew G. Knepley       for (j = loff; j < ldof + loff; ++j) {
4954bbf5ea8SMatthew G. Knepley         for (l = 0; l < bs; ++l) {
4964bbf5ea8SMatthew G. Knepley           PetscInt dof = bs*j + l + subspaceOffset;
4974bbf5ea8SMatthew G. Knepley           PetscHashIAdd(dofs, dof, 0);
4984bbf5ea8SMatthew G. Knepley         }
4994bbf5ea8SMatthew G. Knepley       }
5004bbf5ea8SMatthew G. Knepley     }
5014bbf5ea8SMatthew G. Knepley   }
5024bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
5034bbf5ea8SMatthew G. Knepley }
5044bbf5ea8SMatthew G. Knepley 
5054bbf5ea8SMatthew G. Knepley /* Given two hash tables A and B, compute the keys in B that are not in A, and put them in C */
5064bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchComputeSetDifference_Private(PetscHashI A, PetscHashI B, PetscHashI C)
5074bbf5ea8SMatthew G. Knepley {
5084bbf5ea8SMatthew G. Knepley   PetscHashIIter hi;
5094bbf5ea8SMatthew G. Knepley   PetscInt       key, val;
5104bbf5ea8SMatthew G. Knepley   PetscBool      flg;
5114bbf5ea8SMatthew G. Knepley 
5124bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
5134bbf5ea8SMatthew G. Knepley   PetscHashIClear(C);
5144bbf5ea8SMatthew G. Knepley   PetscHashIIterBegin(B, hi);
5154bbf5ea8SMatthew G. Knepley   while (!PetscHashIIterAtEnd(B, hi)) {
5164bbf5ea8SMatthew G. Knepley     PetscHashIIterGetKeyVal(B, hi, key, val);
5174bbf5ea8SMatthew G. Knepley     PetscHashIIterNext(B, hi);
5184bbf5ea8SMatthew G. Knepley     PetscHashIHasKey(A, key, flg);
5194bbf5ea8SMatthew G. Knepley     if (!flg) {PetscHashIAdd(C, key, val);}
5204bbf5ea8SMatthew G. Knepley   }
5214bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
5224bbf5ea8SMatthew G. Knepley }
5234bbf5ea8SMatthew G. Knepley 
5244bbf5ea8SMatthew G. Knepley /*
5254bbf5ea8SMatthew G. Knepley  * PCPatchCreateCellPatches - create patches.
5264bbf5ea8SMatthew G. Knepley  *
5274bbf5ea8SMatthew G. Knepley  * Input Parameters:
5284bbf5ea8SMatthew G. Knepley  * + dm - The DMPlex object defining the mesh
5294bbf5ea8SMatthew G. Knepley  *
5304bbf5ea8SMatthew G. Knepley  * Output Parameters:
5314bbf5ea8SMatthew G. Knepley  * + cellCounts - Section with counts of cells around each vertex
5324bbf5ea8SMatthew G. Knepley  * - cells - IS of the cell point indices of cells in each patch
5334bbf5ea8SMatthew G. Knepley  */
5344bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchCreateCellPatches(PC pc)
5354bbf5ea8SMatthew G. Knepley {
5364bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
5374bbf5ea8SMatthew G. Knepley   DM             dm, plex;
5384bbf5ea8SMatthew G. Knepley   DMLabel        ghost;
5394bbf5ea8SMatthew G. Knepley   PetscHashI     ht, cht;
5404bbf5ea8SMatthew G. Knepley   PetscInt      *cellsArray = NULL;
5414bbf5ea8SMatthew G. Knepley   PetscInt       numCells;
5424bbf5ea8SMatthew G. Knepley   PetscSection   cellCounts;
5434bbf5ea8SMatthew G. Knepley   PetscInt       pStart, pEnd, vStart, vEnd, v, cStart, cEnd;
5444bbf5ea8SMatthew G. Knepley   PetscBool      flg;
5454bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
5464bbf5ea8SMatthew G. Knepley 
5474bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
5484bbf5ea8SMatthew G. Knepley   /* Used to keep track of the cells in the patch. */
5494bbf5ea8SMatthew G. Knepley   PetscHashICreate(ht);
5504bbf5ea8SMatthew G. Knepley   PetscHashICreate(cht);
5514bbf5ea8SMatthew G. Knepley 
5524bbf5ea8SMatthew G. Knepley   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
5534bbf5ea8SMatthew G. Knepley   if (!dm) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "DM not yet set on patch PC\n");
5544bbf5ea8SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
5554bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetChart(plex, &pStart, &pEnd);CHKERRQ(ierr);
5564bbf5ea8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
5574bbf5ea8SMatthew G. Knepley 
5584bbf5ea8SMatthew G. Knepley   if (patch->user_patches) {
5594bbf5ea8SMatthew G. Knepley     /* compute patch->nuserIS, patch->userIS here */
5604bbf5ea8SMatthew G. Knepley     ierr = patch->userpatchconstructionop(pc, &patch->nuserIS, &patch->userIS, &patch->iterationSet, patch->userpatchconstructctx);CHKERRQ(ierr);
5614bbf5ea8SMatthew G. Knepley     vStart = 0;
5624bbf5ea8SMatthew G. Knepley     vEnd   = patch->nuserIS;
5634bbf5ea8SMatthew G. Knepley   } else if (patch->codim < 0) { /* codim unset */
5644bbf5ea8SMatthew G. Knepley     if (patch->dim < 0) { /* dim unset */
5654bbf5ea8SMatthew G. Knepley       ierr = DMPlexGetDepthStratum(plex, 0, &vStart, &vEnd);CHKERRQ(ierr);
5664bbf5ea8SMatthew G. Knepley     } else { /* dim set */
5674bbf5ea8SMatthew G. Knepley       ierr = DMPlexGetDepthStratum(plex, patch->dim, &vStart, &vEnd);CHKERRQ(ierr);
5684bbf5ea8SMatthew G. Knepley     }
5694bbf5ea8SMatthew G. Knepley   } else { /* codim set */
5704bbf5ea8SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(plex, patch->codim, &vStart, &vEnd);CHKERRQ(ierr);
5714bbf5ea8SMatthew G. Knepley   }
5724bbf5ea8SMatthew G. Knepley 
5734bbf5ea8SMatthew G. Knepley   /* These labels mark the owned points.  We only create patches around points that this process owns. */
5744bbf5ea8SMatthew G. Knepley   /* Replace this with a check of the SF */
5754bbf5ea8SMatthew G. Knepley   ierr = DMGetLabel(dm, "pyop2_ghost", &ghost);CHKERRQ(ierr);
5764bbf5ea8SMatthew G. Knepley   ierr = DMLabelCreateIndex(ghost, pStart, pEnd);CHKERRQ(ierr);
5774bbf5ea8SMatthew G. Knepley 
5784bbf5ea8SMatthew G. Knepley   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->cellCounts);CHKERRQ(ierr);
5794bbf5ea8SMatthew G. Knepley   cellCounts = patch->cellCounts;
5804bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetChart(cellCounts, vStart, vEnd);CHKERRQ(ierr);
5814bbf5ea8SMatthew G. Knepley   /* Count cells in the patch surrounding each entity */
5824bbf5ea8SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
5834bbf5ea8SMatthew G. Knepley     PetscHashIIter hi;
5844bbf5ea8SMatthew G. Knepley     PetscInt       chtSize;
5854bbf5ea8SMatthew G. Knepley 
5864bbf5ea8SMatthew G. Knepley     if (!patch->user_patches) {
5874bbf5ea8SMatthew G. Knepley       ierr = DMLabelHasPoint(ghost, v, &flg);CHKERRQ(ierr);
5884bbf5ea8SMatthew G. Knepley       /* Not an owned entity, don't make a cell patch. */
5894bbf5ea8SMatthew G. Knepley       if (flg) continue;
5904bbf5ea8SMatthew G. Knepley     }
5914bbf5ea8SMatthew G. Knepley 
5924bbf5ea8SMatthew G. Knepley     ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr);
5934bbf5ea8SMatthew G. Knepley     ierr = PCPatchCompleteCellPatch(dm, ht, cht);CHKERRQ(ierr);
5944bbf5ea8SMatthew G. Knepley     PetscHashISize(cht, chtSize);
5954bbf5ea8SMatthew G. Knepley     /* empty patch, continue */
5964bbf5ea8SMatthew G. Knepley     if (chtSize == 0) continue;
5974bbf5ea8SMatthew G. Knepley 
5984bbf5ea8SMatthew G. Knepley     /* safe because size(cht) > 0 from above */
5994bbf5ea8SMatthew G. Knepley     PetscHashIIterBegin(cht, hi);
6004bbf5ea8SMatthew G. Knepley     while (!PetscHashIIterAtEnd(cht, hi)) {
6014bbf5ea8SMatthew G. Knepley       PetscInt point;
6024bbf5ea8SMatthew G. Knepley 
6034bbf5ea8SMatthew G. Knepley       PetscHashIIterGetKey(cht, hi, point);
6044bbf5ea8SMatthew G. Knepley       if (point >= cStart && point < cEnd) {
6054bbf5ea8SMatthew G. Knepley         ierr = PetscSectionAddDof(cellCounts, v, 1);CHKERRQ(ierr);
6064bbf5ea8SMatthew G. Knepley       }
6074bbf5ea8SMatthew G. Knepley       PetscHashIIterNext(cht, hi);
6084bbf5ea8SMatthew G. Knepley     }
6094bbf5ea8SMatthew G. Knepley   }
6104bbf5ea8SMatthew G. Knepley   ierr = DMLabelDestroyIndex(ghost);CHKERRQ(ierr);
6114bbf5ea8SMatthew G. Knepley 
6124bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetUp(cellCounts);CHKERRQ(ierr);
6134bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr);
6144bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numCells, &cellsArray);CHKERRQ(ierr);
6154bbf5ea8SMatthew G. Knepley 
6164bbf5ea8SMatthew G. Knepley   /* Now that we know how much space we need, run through again and actually remember the cells. */
6174bbf5ea8SMatthew G. Knepley   for (v = vStart; v < vEnd; v++ ) {
6184bbf5ea8SMatthew G. Knepley     PetscHashIIter hi;
6194bbf5ea8SMatthew G. Knepley     PetscInt       ndof, off;
6204bbf5ea8SMatthew G. Knepley 
6214bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetDof(cellCounts, v, &ndof);CHKERRQ(ierr);
6224bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
6234bbf5ea8SMatthew G. Knepley     if (ndof <= 0) continue;
6244bbf5ea8SMatthew G. Knepley     ierr = patch->patchconstructop((void *) patch, dm, v, ht);CHKERRQ(ierr);
6254bbf5ea8SMatthew G. Knepley     ierr = PCPatchCompleteCellPatch(dm, ht, cht);CHKERRQ(ierr);
6264bbf5ea8SMatthew G. Knepley     ndof = 0;
6274bbf5ea8SMatthew G. Knepley     PetscHashIIterBegin(cht, hi);
6284bbf5ea8SMatthew G. Knepley     while (!PetscHashIIterAtEnd(cht, hi)) {
6294bbf5ea8SMatthew G. Knepley       PetscInt point;
6304bbf5ea8SMatthew G. Knepley 
6314bbf5ea8SMatthew G. Knepley       PetscHashIIterGetKey(cht, hi, point);
6324bbf5ea8SMatthew G. Knepley       if (point >= cStart && point < cEnd) {
6334bbf5ea8SMatthew G. Knepley         cellsArray[ndof + off] = point;
6344bbf5ea8SMatthew G. Knepley         ++ndof;
6354bbf5ea8SMatthew G. Knepley       }
6364bbf5ea8SMatthew G. Knepley       PetscHashIIterNext(cht, hi);
6374bbf5ea8SMatthew G. Knepley     }
6384bbf5ea8SMatthew G. Knepley   }
6394bbf5ea8SMatthew G. Knepley 
6404bbf5ea8SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numCells, cellsArray, PETSC_OWN_POINTER, &patch->cells);CHKERRQ(ierr);
6414bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
6424bbf5ea8SMatthew G. Knepley   patch->npatch = pEnd - pStart;
6434bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(ht);
6444bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(cht);
6454bbf5ea8SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
6464bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
6474bbf5ea8SMatthew G. Knepley }
6484bbf5ea8SMatthew G. Knepley 
6494bbf5ea8SMatthew G. Knepley /*
6504bbf5ea8SMatthew G. Knepley  * PCPatchCreateCellPatchDiscretisationInfo - Build the dof maps for cell patches
6514bbf5ea8SMatthew G. Knepley  *
6524bbf5ea8SMatthew G. Knepley  * Input Parameters:
6534bbf5ea8SMatthew G. Knepley  * + dm - The DMPlex object defining the mesh
6544bbf5ea8SMatthew G. Knepley  * . cellCounts - Section with counts of cells around each vertex
6554bbf5ea8SMatthew G. Knepley  * . cells - IS of the cell point indices of cells in each patch
6564bbf5ea8SMatthew G. Knepley  * . cellNumbering - Section mapping plex cell points to Firedrake cell indices.
6574bbf5ea8SMatthew G. Knepley  * . nodesPerCell - number of nodes per cell.
6584bbf5ea8SMatthew G. Knepley  * - cellNodeMap - map from cells to node indices (nodesPerCell * numCells)
6594bbf5ea8SMatthew G. Knepley  *
6604bbf5ea8SMatthew G. Knepley  * Output Parameters:
6614bbf5ea8SMatthew G. Knepley  * + dofs - IS of local dof numbers of each cell in the patch
6624bbf5ea8SMatthew G. Knepley  * . gtolCounts - Section with counts of dofs per cell patch
6634bbf5ea8SMatthew G. Knepley  * - gtol - IS mapping from global dofs to local dofs for each patch.
6644bbf5ea8SMatthew G. Knepley  */
6654bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchCreateCellPatchDiscretisationInfo(PC pc)
6664bbf5ea8SMatthew G. Knepley {
6674bbf5ea8SMatthew G. Knepley   PC_PATCH       *patch           = (PC_PATCH *) pc->data;
6684bbf5ea8SMatthew G. Knepley   PetscSection    cellCounts      = patch->cellCounts;
6694bbf5ea8SMatthew G. Knepley   PetscSection    gtolCounts;
6704bbf5ea8SMatthew G. Knepley   IS              cells           = patch->cells;
6714bbf5ea8SMatthew G. Knepley   PetscSection    cellNumbering   = patch->cellNumbering;
6724bbf5ea8SMatthew G. Knepley   PetscInt        numCells;
6734bbf5ea8SMatthew G. Knepley   PetscInt        numDofs;
6744bbf5ea8SMatthew G. Knepley   PetscInt        numGlobalDofs;
6754bbf5ea8SMatthew G. Knepley   PetscInt        totalDofsPerCell = patch->totalDofsPerCell;
6764bbf5ea8SMatthew G. Knepley   PetscInt        vStart, vEnd, v;
6774bbf5ea8SMatthew G. Knepley   const PetscInt *cellsArray;
6784bbf5ea8SMatthew G. Knepley   PetscInt       *newCellsArray   = NULL;
6794bbf5ea8SMatthew G. Knepley   PetscInt       *dofsArray       = NULL;
6804bbf5ea8SMatthew G. Knepley   PetscInt       *asmArray        = NULL;
6814bbf5ea8SMatthew G. Knepley   PetscInt       *globalDofsArray = NULL;
6824bbf5ea8SMatthew G. Knepley   PetscInt        globalIndex     = 0;
6834bbf5ea8SMatthew G. Knepley   PetscInt        key             = 0;
6844bbf5ea8SMatthew G. Knepley   PetscInt        asmKey          = 0;
6854bbf5ea8SMatthew G. Knepley   PetscHashI      ht;
6864bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
6874bbf5ea8SMatthew G. Knepley 
6884bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
6894bbf5ea8SMatthew G. Knepley   /* dofcounts section is cellcounts section * dofPerCell */
6904bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(cellCounts, &numCells);CHKERRQ(ierr);
6914bbf5ea8SMatthew G. Knepley   numDofs = numCells * totalDofsPerCell;
6924bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numDofs, &dofsArray);CHKERRQ(ierr);
6934bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numDofs, &asmArray);CHKERRQ(ierr);
6944bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numCells, &newCellsArray);CHKERRQ(ierr);
6954bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetChart(cellCounts, &vStart, &vEnd);CHKERRQ(ierr);
6964bbf5ea8SMatthew G. Knepley   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->gtolCounts);CHKERRQ(ierr);
6974bbf5ea8SMatthew G. Knepley   gtolCounts = patch->gtolCounts;
6984bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetChart(gtolCounts, vStart, vEnd);CHKERRQ(ierr);
6994bbf5ea8SMatthew G. Knepley 
7004bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(cells, &cellsArray);CHKERRQ(ierr);
7014bbf5ea8SMatthew G. Knepley   PetscHashICreate(ht);
7024bbf5ea8SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
7034bbf5ea8SMatthew G. Knepley     PetscInt localIndex = 0;
7044bbf5ea8SMatthew G. Knepley     PetscInt dof, off, i, j, k, l;
7054bbf5ea8SMatthew G. Knepley 
7064bbf5ea8SMatthew G. Knepley     PetscHashIClear(ht);
7074bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
7084bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
7094bbf5ea8SMatthew G. Knepley     if (dof <= 0) continue;
7104bbf5ea8SMatthew G. Knepley 
7114bbf5ea8SMatthew G. Knepley     for (k = 0; k < patch->nsubspaces; ++k) {
7124bbf5ea8SMatthew G. Knepley       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
7134bbf5ea8SMatthew G. Knepley       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
7144bbf5ea8SMatthew G. Knepley       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
7154bbf5ea8SMatthew G. Knepley       PetscInt        bs             = patch->bs[k];
7164bbf5ea8SMatthew G. Knepley 
7174bbf5ea8SMatthew G. Knepley       for (i = off; i < off + dof; ++i) {
7184bbf5ea8SMatthew G. Knepley         /* Walk over the cells in this patch. */
7194bbf5ea8SMatthew G. Knepley         const PetscInt c = cellsArray[i];
7204bbf5ea8SMatthew G. Knepley         PetscInt       cell;
7214bbf5ea8SMatthew G. Knepley 
7224bbf5ea8SMatthew G. Knepley         ierr = PetscSectionGetDof(cellNumbering, c, &cell);CHKERRQ(ierr);
7234bbf5ea8SMatthew G. Knepley         if (cell <= 0) SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_OUTOFRANGE, "Cell %D doesn't appear in cell numbering map", c);
7244bbf5ea8SMatthew G. Knepley         ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);
7254bbf5ea8SMatthew G. Knepley         newCellsArray[i] = cell;
7264bbf5ea8SMatthew G. Knepley         for (j = 0; j < nodesPerCell; ++j) {
7274bbf5ea8SMatthew G. Knepley           /* For each global dof, map it into contiguous local storage. */
7284bbf5ea8SMatthew G. Knepley           const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset;
7294bbf5ea8SMatthew G. Knepley           /* finally, loop over block size */
7304bbf5ea8SMatthew G. Knepley           for (l = 0; l < bs; ++l) {
7314bbf5ea8SMatthew G. Knepley             PetscInt localDof;
7324bbf5ea8SMatthew G. Knepley 
7334bbf5ea8SMatthew G. Knepley             PetscHashIMap(ht, globalDof + l, localDof);
7344bbf5ea8SMatthew G. Knepley             if (localDof == -1) {
7354bbf5ea8SMatthew G. Knepley               localDof = localIndex++;
7364bbf5ea8SMatthew G. Knepley               PetscHashIAdd(ht, globalDof + l, localDof);
7374bbf5ea8SMatthew G. Knepley             }
7384bbf5ea8SMatthew G. Knepley             if (globalIndex >= numDofs) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Found more dofs %D than expected %D", globalIndex+1, numDofs);
7394bbf5ea8SMatthew G. Knepley             /* And store. */
7404bbf5ea8SMatthew G. Knepley             dofsArray[globalIndex++] = localDof;
7414bbf5ea8SMatthew G. Knepley           }
7424bbf5ea8SMatthew G. Knepley         }
7434bbf5ea8SMatthew G. Knepley       }
7444bbf5ea8SMatthew G. Knepley     }
7454bbf5ea8SMatthew G. Knepley     /* How many local dofs in this patch? */
7464bbf5ea8SMatthew G. Knepley     PetscHashISize(ht, dof);
7474bbf5ea8SMatthew G. Knepley     ierr = PetscSectionSetDof(gtolCounts, v, dof);CHKERRQ(ierr);
7484bbf5ea8SMatthew G. Knepley   }
7494bbf5ea8SMatthew G. Knepley   if (globalIndex != numDofs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Expected number of dofs (%d) doesn't match found number (%d)", numDofs, globalIndex);
7504bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetUp(gtolCounts);CHKERRQ(ierr);
7514bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(gtolCounts, &numGlobalDofs);CHKERRQ(ierr);
7524bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(numGlobalDofs, &globalDofsArray);CHKERRQ(ierr);
7534bbf5ea8SMatthew G. Knepley 
7544bbf5ea8SMatthew G. Knepley   /* Now populate the global to local map.  This could be merged into the above loop if we were willing to deal with reallocs. */
7554bbf5ea8SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
7564bbf5ea8SMatthew G. Knepley     PetscHashIIter hi;
7574bbf5ea8SMatthew G. Knepley     PetscInt       dof, off, i, j, k, l;
7584bbf5ea8SMatthew G. Knepley 
7594bbf5ea8SMatthew G. Knepley     PetscHashIClear(ht);
7604bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetDof(cellCounts, v, &dof);CHKERRQ(ierr);
7614bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetOffset(cellCounts, v, &off);CHKERRQ(ierr);
7624bbf5ea8SMatthew G. Knepley     if (dof <= 0) continue;
7634bbf5ea8SMatthew G. Knepley 
7644bbf5ea8SMatthew G. Knepley     for (k = 0; k < patch->nsubspaces; ++k) {
7654bbf5ea8SMatthew G. Knepley       const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
7664bbf5ea8SMatthew G. Knepley       PetscInt        nodesPerCell   = patch->nodesPerCell[k];
7674bbf5ea8SMatthew G. Knepley       PetscInt        subspaceOffset = patch->subspaceOffsets[k];
7684bbf5ea8SMatthew G. Knepley       PetscInt        bs             = patch->bs[k];
7694bbf5ea8SMatthew G. Knepley 
7704bbf5ea8SMatthew G. Knepley       for (i = off; i < off + dof; ++i) {
7714bbf5ea8SMatthew G. Knepley         /* Reconstruct mapping of global-to-local on this patch. */
7724bbf5ea8SMatthew G. Knepley         const PetscInt c = cellsArray[i];
7734bbf5ea8SMatthew G. Knepley         PetscInt       cell;
7744bbf5ea8SMatthew G. Knepley 
7754bbf5ea8SMatthew G. Knepley         ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);
7764bbf5ea8SMatthew G. Knepley         for (j = 0; j < nodesPerCell; ++j) {
7774bbf5ea8SMatthew G. Knepley           for (l = 0; l < bs; ++l) {
7784bbf5ea8SMatthew G. Knepley             const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset + l;
7794bbf5ea8SMatthew G. Knepley             const PetscInt localDof  = dofsArray[key];
7804bbf5ea8SMatthew G. Knepley 
7814bbf5ea8SMatthew G. Knepley             key += 1;
7824bbf5ea8SMatthew G. Knepley             PetscHashIAdd(ht, globalDof, localDof);
7834bbf5ea8SMatthew G. Knepley           }
7844bbf5ea8SMatthew G. Knepley         }
7854bbf5ea8SMatthew G. Knepley       }
7864bbf5ea8SMatthew G. Knepley       if (dof > 0) {
7874bbf5ea8SMatthew G. Knepley         /* Shove it in the output data structure. */
7884bbf5ea8SMatthew G. Knepley         PetscInt goff;
7894bbf5ea8SMatthew G. Knepley 
7904bbf5ea8SMatthew G. Knepley         ierr = PetscSectionGetOffset(gtolCounts, v, &goff);CHKERRQ(ierr);
7914bbf5ea8SMatthew G. Knepley         PetscHashIIterBegin(ht, hi);
7924bbf5ea8SMatthew G. Knepley         while (!PetscHashIIterAtEnd(ht, hi)) {
7934bbf5ea8SMatthew G. Knepley           PetscInt globalDof, localDof;
7944bbf5ea8SMatthew G. Knepley 
7954bbf5ea8SMatthew G. Knepley           PetscHashIIterGetKeyVal(ht, hi, globalDof, localDof);
7964bbf5ea8SMatthew G. Knepley           if (globalDof >= 0) globalDofsArray[goff + localDof] = globalDof;
7974bbf5ea8SMatthew G. Knepley           PetscHashIIterNext(ht, hi);
7984bbf5ea8SMatthew G. Knepley         }
7994bbf5ea8SMatthew G. Knepley       }
8004bbf5ea8SMatthew G. Knepley     }
8014bbf5ea8SMatthew G. Knepley 
8024bbf5ea8SMatthew G. Knepley     /* At this point, we have a hash table ht built that maps globalDof -> localDof.
8034bbf5ea8SMatthew G. Knepley      We need to create the dof table laid out cellwise first, then by subspace,
8044bbf5ea8SMatthew G. Knepley      as the assembler assembles cell-wise and we need to stuff the different
8054bbf5ea8SMatthew G. Knepley      contributions of the different function spaces to the right places. So we loop
8064bbf5ea8SMatthew G. Knepley      over cells, then over subspaces. */
8074bbf5ea8SMatthew G. Knepley     if (patch->nsubspaces > 1) { /* for nsubspaces = 1, data we need is already in dofsArray */
8084bbf5ea8SMatthew G. Knepley       for (i = off; i < off + dof; ++i) {
8094bbf5ea8SMatthew G. Knepley         const PetscInt c = cellsArray[i];
8104bbf5ea8SMatthew G. Knepley         PetscInt       cell;
8114bbf5ea8SMatthew G. Knepley 
8124bbf5ea8SMatthew G. Knepley         ierr = PetscSectionGetOffset(cellNumbering, c, &cell);CHKERRQ(ierr);
8134bbf5ea8SMatthew G. Knepley         for (k = 0; k < patch->nsubspaces; ++k) {
8144bbf5ea8SMatthew G. Knepley           const PetscInt *cellNodeMap    = patch->cellNodeMap[k];
8154bbf5ea8SMatthew G. Knepley           PetscInt        nodesPerCell   = patch->nodesPerCell[k];
8164bbf5ea8SMatthew G. Knepley           PetscInt        subspaceOffset = patch->subspaceOffsets[k];
8174bbf5ea8SMatthew G. Knepley           PetscInt        bs             = patch->bs[k];
8184bbf5ea8SMatthew G. Knepley 
8194bbf5ea8SMatthew G. Knepley           for (j = 0; j < nodesPerCell; ++j) {
8204bbf5ea8SMatthew G. Knepley             for (l = 0; l < bs; ++l) {
8214bbf5ea8SMatthew G. Knepley               const PetscInt globalDof = cellNodeMap[cell*nodesPerCell + j]*bs + subspaceOffset + l;
8224bbf5ea8SMatthew G. Knepley               PetscInt       localDof;
8234bbf5ea8SMatthew G. Knepley 
8244bbf5ea8SMatthew G. Knepley               PetscHashIMap(ht, globalDof, localDof);
8254bbf5ea8SMatthew G. Knepley               asmArray[asmKey++] = localDof;
8264bbf5ea8SMatthew G. Knepley             }
8274bbf5ea8SMatthew G. Knepley           }
8284bbf5ea8SMatthew G. Knepley         }
8294bbf5ea8SMatthew G. Knepley       }
8304bbf5ea8SMatthew G. Knepley     }
8314bbf5ea8SMatthew G. Knepley   }
8324bbf5ea8SMatthew G. Knepley   if (1 == patch->nsubspaces) {ierr = PetscMemcpy(asmArray, dofsArray, numDofs * sizeof(PetscInt));CHKERRQ(ierr);}
8334bbf5ea8SMatthew G. Knepley 
8344bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(ht);
8354bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(cells, &cellsArray);CHKERRQ(ierr);
8364bbf5ea8SMatthew G. Knepley   ierr = PetscFree(dofsArray);CHKERRQ(ierr);
8374bbf5ea8SMatthew G. Knepley   /* Replace cell indices with firedrake-numbered ones. */
8384bbf5ea8SMatthew G. Knepley   ierr = ISGeneralSetIndices(cells, numCells, (const PetscInt *) newCellsArray, PETSC_OWN_POINTER);CHKERRQ(ierr);
8394bbf5ea8SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numGlobalDofs, globalDofsArray, PETSC_OWN_POINTER, &patch->gtol);CHKERRQ(ierr);
8404bbf5ea8SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, numDofs, asmArray, PETSC_OWN_POINTER, &patch->dofs);CHKERRQ(ierr);
8414bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
8424bbf5ea8SMatthew G. Knepley }
8434bbf5ea8SMatthew G. Knepley 
8444bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchCreateCellPatchBCs(PC pc)
8454bbf5ea8SMatthew G. Knepley {
8464bbf5ea8SMatthew G. Knepley   PC_PATCH       *patch      = (PC_PATCH *) pc->data;
8474bbf5ea8SMatthew G. Knepley   const PetscInt *bcNodes    = NULL;
8484bbf5ea8SMatthew G. Knepley   PetscSection    gtolCounts = patch->gtolCounts;
8494bbf5ea8SMatthew G. Knepley   IS              gtol       = patch->gtol;
8504bbf5ea8SMatthew G. Knepley   DM              dm;
8514bbf5ea8SMatthew G. Knepley   PetscInt        numBcs;
8524bbf5ea8SMatthew G. Knepley   PetscSection    bcCounts;
8534bbf5ea8SMatthew G. Knepley   PetscHashI      globalBcs, localBcs, patchDofs;
8544bbf5ea8SMatthew G. Knepley   PetscHashI      ownedpts, seenpts, owneddofs, seendofs, artificialbcs;
8554bbf5ea8SMatthew G. Knepley   PetscHashIIter  hi;
8564bbf5ea8SMatthew G. Knepley   PetscInt       *bcsArray     = NULL;
8574bbf5ea8SMatthew G. Knepley   PetscInt       *multBcsArray = NULL;
8584bbf5ea8SMatthew G. Knepley   const PetscInt *gtolArray;
8594bbf5ea8SMatthew G. Knepley   PetscInt        vStart, vEnd, v, i;
8604bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
8614bbf5ea8SMatthew G. Knepley 
8624bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
8634bbf5ea8SMatthew G. Knepley   ierr = PCGetDM(pc, &dm);CHKERRQ(ierr);
8644bbf5ea8SMatthew G. Knepley   PetscHashICreate(globalBcs);
8654bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(patch->ghostBcNodes, &bcNodes);CHKERRQ(ierr);
8664bbf5ea8SMatthew G. Knepley   ierr = ISGetSize(patch->ghostBcNodes, &numBcs);CHKERRQ(ierr);
8674bbf5ea8SMatthew G. Knepley   for (i = 0; i < numBcs; ++i) {
8684bbf5ea8SMatthew G. Knepley     /* these are already in concatenated numbering */
8694bbf5ea8SMatthew G. Knepley     PetscHashIAdd(globalBcs, bcNodes[i], 0);
8704bbf5ea8SMatthew G. Knepley   }
8714bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(patch->ghostBcNodes, &bcNodes);CHKERRQ(ierr);
8724bbf5ea8SMatthew G. Knepley   PetscHashICreate(patchDofs);
8734bbf5ea8SMatthew G. Knepley   PetscHashICreate(localBcs);
8744bbf5ea8SMatthew G. Knepley   PetscHashICreate(ownedpts);
8754bbf5ea8SMatthew G. Knepley   PetscHashICreate(seenpts);
8764bbf5ea8SMatthew G. Knepley   PetscHashICreate(owneddofs);
8774bbf5ea8SMatthew G. Knepley   PetscHashICreate(seendofs);
8784bbf5ea8SMatthew G. Knepley   PetscHashICreate(artificialbcs);
8794bbf5ea8SMatthew G. Knepley 
8804bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetChart(patch->cellCounts, &vStart, &vEnd);CHKERRQ(ierr);
8814bbf5ea8SMatthew G. Knepley   ierr = PetscSectionCreate(PETSC_COMM_SELF, &patch->bcCounts);CHKERRQ(ierr);
8824bbf5ea8SMatthew G. Knepley   bcCounts = patch->bcCounts;
8834bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetChart(bcCounts, vStart, vEnd);CHKERRQ(ierr);
8844bbf5ea8SMatthew G. Knepley   ierr = PetscMalloc1(vEnd - vStart, &patch->bcs);CHKERRQ(ierr);
8854bbf5ea8SMatthew G. Knepley   if (patch->multiplicative) {ierr = PetscMalloc1(vEnd - vStart, &patch->multBcs);CHKERRQ(ierr);}
8864bbf5ea8SMatthew G. Knepley 
8874bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(gtol, &gtolArray);CHKERRQ(ierr);
8884bbf5ea8SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
8894bbf5ea8SMatthew G. Knepley     PetscInt bcIndex     = 0;
8904bbf5ea8SMatthew G. Knepley     PetscInt multBcIndex = 0;
8914bbf5ea8SMatthew G. Knepley     PetscInt numBcs, dof, off;
8924bbf5ea8SMatthew G. Knepley 
8934bbf5ea8SMatthew G. Knepley     PetscHashIClear(patchDofs);
8944bbf5ea8SMatthew G. Knepley     PetscHashIClear(localBcs);
8954bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetDof(gtolCounts, v, &dof);CHKERRQ(ierr);
8964bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetOffset(gtolCounts, v, &off);CHKERRQ(ierr);
8974bbf5ea8SMatthew G. Knepley 
8984bbf5ea8SMatthew G. Knepley     if (dof <= 0) {
8994bbf5ea8SMatthew G. Knepley       patch->bcs[v-vStart] = NULL;
9004bbf5ea8SMatthew G. Knepley       if (patch->multiplicative) patch->multBcs[v-vStart] = NULL;
9014bbf5ea8SMatthew G. Knepley       continue;
9024bbf5ea8SMatthew G. Knepley     }
9034bbf5ea8SMatthew G. Knepley 
9044bbf5ea8SMatthew G. Knepley     for (i = off; i < off + dof; ++i) {
9054bbf5ea8SMatthew G. Knepley       const PetscInt globalDof = gtolArray[i];
9064bbf5ea8SMatthew G. Knepley       const PetscInt localDof  = i-off;
9074bbf5ea8SMatthew G. Knepley       PetscBool      flg;
9084bbf5ea8SMatthew G. Knepley 
9094bbf5ea8SMatthew G. Knepley       PetscHashIAdd(patchDofs, globalDof, localDof);
9104bbf5ea8SMatthew G. Knepley       PetscHashIHasKey(globalBcs, globalDof, flg);
9114bbf5ea8SMatthew G. Knepley       if (flg) {PetscHashIAdd(localBcs, localDof, 0);}
9124bbf5ea8SMatthew G. Knepley     }
9134bbf5ea8SMatthew G. Knepley 
9144bbf5ea8SMatthew G. Knepley     /* If we're doing multiplicative, make the BC data structures now corresponding solely to actual globally imposed Dirichlet BCs */
9154bbf5ea8SMatthew G. Knepley     if (patch->multiplicative) {
9164bbf5ea8SMatthew G. Knepley       PetscHashISize(localBcs, numBcs);
9174bbf5ea8SMatthew G. Knepley       ierr = PetscMalloc1(numBcs, &multBcsArray);CHKERRQ(ierr);
9184bbf5ea8SMatthew G. Knepley       PetscHashIGetKeys(localBcs, &multBcIndex, multBcsArray);
9194bbf5ea8SMatthew G. Knepley       ierr = PetscSortInt(numBcs, multBcsArray);CHKERRQ(ierr);
9204bbf5ea8SMatthew G. Knepley       ierr = ISCreateGeneral(PETSC_COMM_SELF, numBcs, multBcsArray, PETSC_OWN_POINTER, &(patch->multBcs[v-vStart]));CHKERRQ(ierr);
9214bbf5ea8SMatthew G. Knepley     }
9224bbf5ea8SMatthew G. Knepley 
9234bbf5ea8SMatthew G. Knepley     /* Now figure out the artificial BCs: the set difference of {dofs on entities I see on the patch}\{dofs I am responsible for updating} */
9244bbf5ea8SMatthew G. Knepley     ierr = patch->patchconstructop((void*) patch, dm, v, ownedpts);CHKERRQ(ierr);
9254bbf5ea8SMatthew G. Knepley     ierr = PCPatchCompleteCellPatch(dm, ownedpts, seenpts);CHKERRQ(ierr);
9264bbf5ea8SMatthew G. Knepley     ierr = PCPatchGetPointDofs(patch, ownedpts, owneddofs, v, patch->exclude_subspace);CHKERRQ(ierr);
9274bbf5ea8SMatthew G. Knepley     ierr = PCPatchGetPointDofs(patch, seenpts, seendofs, v, -1);CHKERRQ(ierr);
9284bbf5ea8SMatthew G. Knepley     ierr = PCPatchComputeSetDifference_Private(owneddofs, seendofs, artificialbcs);CHKERRQ(ierr);
9294bbf5ea8SMatthew G. Knepley 
9304bbf5ea8SMatthew G. Knepley     if (patch->print_patches) {
9314bbf5ea8SMatthew G. Knepley       PetscHashI globalbcdofs;
9324bbf5ea8SMatthew G. Knepley       MPI_Comm   comm;
9334bbf5ea8SMatthew G. Knepley 
9344bbf5ea8SMatthew G. Knepley       PetscHashICreate(globalbcdofs);
9354bbf5ea8SMatthew G. Knepley 
9364bbf5ea8SMatthew G. Knepley       ierr = PetscObjectGetComm((PetscObject) pc, &comm);CHKERRQ(ierr);
9374bbf5ea8SMatthew G. Knepley       ierr = PetscSynchronizedPrintf(comm, "Patch %d: owned dofs:\n", v);CHKERRQ(ierr);
9384bbf5ea8SMatthew G. Knepley       PetscHashIIterBegin(owneddofs, hi);
9394bbf5ea8SMatthew G. Knepley       while (!PetscHashIIterAtEnd(owneddofs, hi)) {
9404bbf5ea8SMatthew G. Knepley         PetscInt globalDof;
9414bbf5ea8SMatthew G. Knepley 
9424bbf5ea8SMatthew G. Knepley         PetscHashIIterGetKey(owneddofs, hi, globalDof);
9434bbf5ea8SMatthew G. Knepley         PetscHashIIterNext(owneddofs, hi);
9444bbf5ea8SMatthew G. Knepley         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr);
9454bbf5ea8SMatthew G. Knepley       }
9464bbf5ea8SMatthew G. Knepley       ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr);
9474bbf5ea8SMatthew G. Knepley       ierr = PetscSynchronizedPrintf(comm, "Patch %d: seen dofs:\n", v);CHKERRQ(ierr);
9484bbf5ea8SMatthew G. Knepley       PetscHashIIterBegin(seendofs, hi);
9494bbf5ea8SMatthew G. Knepley       while (!PetscHashIIterAtEnd(seendofs, hi)) {
9504bbf5ea8SMatthew G. Knepley         PetscInt globalDof;
9514bbf5ea8SMatthew G. Knepley         PetscBool flg;
9524bbf5ea8SMatthew G. Knepley 
9534bbf5ea8SMatthew G. Knepley         PetscHashIIterGetKey(seendofs, hi, globalDof);
9544bbf5ea8SMatthew G. Knepley         PetscHashIIterNext(seendofs, hi);
9554bbf5ea8SMatthew G. Knepley         ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr);
9564bbf5ea8SMatthew G. Knepley         PetscHashIHasKey(globalBcs, globalDof, flg);
9574bbf5ea8SMatthew G. Knepley         if (flg) {PetscHashIAdd(globalbcdofs, globalDof, 0);}
9584bbf5ea8SMatthew G. Knepley       }
9594bbf5ea8SMatthew G. Knepley       ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr);
9604bbf5ea8SMatthew G. Knepley       ierr = PetscSynchronizedPrintf(comm, "Patch %d: global BCs:\n", v);CHKERRQ(ierr);
9614bbf5ea8SMatthew G. Knepley       PetscHashISize(globalbcdofs, numBcs);
9624bbf5ea8SMatthew G. Knepley       if (numBcs > 0) {
9634bbf5ea8SMatthew G. Knepley         PetscHashIIterBegin(globalbcdofs, hi);
9644bbf5ea8SMatthew G. Knepley         while (!PetscHashIIterAtEnd(globalbcdofs, hi)) {
9654bbf5ea8SMatthew G. Knepley           PetscInt globalDof;
9664bbf5ea8SMatthew G. Knepley           PetscHashIIterGetKey(globalbcdofs, hi, globalDof);
9674bbf5ea8SMatthew G. Knepley           PetscHashIIterNext(globalbcdofs, hi);
9684bbf5ea8SMatthew G. Knepley           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr);
9694bbf5ea8SMatthew G. Knepley         }
9704bbf5ea8SMatthew G. Knepley       }
9714bbf5ea8SMatthew G. Knepley       ierr = PetscSynchronizedPrintf(comm, "\n");CHKERRQ(ierr);
9724bbf5ea8SMatthew G. Knepley       ierr = PetscSynchronizedPrintf(comm, "Patch %d: artificial BCs:\n", v);CHKERRQ(ierr);
9734bbf5ea8SMatthew G. Knepley       PetscHashISize(artificialbcs, numBcs);
9744bbf5ea8SMatthew G. Knepley       if (numBcs > 0) {
9754bbf5ea8SMatthew G. Knepley         PetscHashIIterBegin(artificialbcs, hi);
9764bbf5ea8SMatthew G. Knepley         while (!PetscHashIIterAtEnd(artificialbcs, hi)) {
9774bbf5ea8SMatthew G. Knepley           PetscInt globalDof;
9784bbf5ea8SMatthew G. Knepley           PetscHashIIterGetKey(artificialbcs, hi, globalDof);
9794bbf5ea8SMatthew G. Knepley           PetscHashIIterNext(artificialbcs, hi);
9804bbf5ea8SMatthew G. Knepley           ierr = PetscSynchronizedPrintf(comm, "%d ", globalDof);CHKERRQ(ierr);
9814bbf5ea8SMatthew G. Knepley         }
9824bbf5ea8SMatthew G. Knepley       }
9834bbf5ea8SMatthew G. Knepley       ierr = PetscSynchronizedPrintf(comm, "\n\n");CHKERRQ(ierr);
9844bbf5ea8SMatthew G. Knepley       ierr = PetscSynchronizedFlush(comm, PETSC_STDOUT);CHKERRQ(ierr);
9854bbf5ea8SMatthew G. Knepley       PetscHashIDestroy(globalbcdofs);
9864bbf5ea8SMatthew G. Knepley     }
9874bbf5ea8SMatthew G. Knepley 
9884bbf5ea8SMatthew G. Knepley     PetscHashISize(artificialbcs, numBcs);
9894bbf5ea8SMatthew G. Knepley     if (numBcs > 0) {
9904bbf5ea8SMatthew G. Knepley       PetscHashIIterBegin(artificialbcs, hi);
9914bbf5ea8SMatthew G. Knepley       while (!PetscHashIIterAtEnd(artificialbcs, hi)) {
9924bbf5ea8SMatthew G. Knepley         PetscInt globalDof, localDof;
9934bbf5ea8SMatthew G. Knepley         PetscHashIIterGetKey(artificialbcs, hi, globalDof);
9944bbf5ea8SMatthew G. Knepley         PetscHashIIterNext(artificialbcs, hi);
9954bbf5ea8SMatthew G. Knepley         PetscHashIMap(patchDofs, globalDof, localDof);
9964bbf5ea8SMatthew G. Knepley         if (localDof == -1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Patch %d Didn't find dof %d in patch\n", v-vStart, globalDof);
9974bbf5ea8SMatthew G. Knepley         PetscHashIAdd(localBcs, localDof, 0);
9984bbf5ea8SMatthew G. Knepley       }
9994bbf5ea8SMatthew G. Knepley     }
10004bbf5ea8SMatthew G. Knepley 
10014bbf5ea8SMatthew G. Knepley     /* OK, now we have a hash table with all the bcs indicated by the artificial and global bcs */
10024bbf5ea8SMatthew G. Knepley     PetscHashISize(localBcs, numBcs);
10034bbf5ea8SMatthew G. Knepley     ierr = PetscSectionSetDof(bcCounts, v, numBcs);CHKERRQ(ierr);
10044bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(numBcs, &bcsArray);CHKERRQ(ierr);
10054bbf5ea8SMatthew G. Knepley     PetscHashIGetKeys(localBcs, &bcIndex, bcsArray);
10064bbf5ea8SMatthew G. Knepley     ierr = PetscSortInt(numBcs, bcsArray);CHKERRQ(ierr);
10074bbf5ea8SMatthew G. Knepley     ierr = ISCreateGeneral(PETSC_COMM_SELF, numBcs, bcsArray, PETSC_OWN_POINTER, &(patch->bcs[v - vStart]));CHKERRQ(ierr);
10084bbf5ea8SMatthew G. Knepley   }
10094bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(gtol, &gtolArray);CHKERRQ(ierr);
10104bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(artificialbcs);
10114bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(seendofs);
10124bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(owneddofs);
10134bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(seenpts);
10144bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(ownedpts);
10154bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(localBcs);
10164bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(patchDofs);
10174bbf5ea8SMatthew G. Knepley   PetscHashIDestroy(globalBcs);
10184bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr);
10194bbf5ea8SMatthew G. Knepley   ierr = PetscSectionSetUp(bcCounts);CHKERRQ(ierr);
10204bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
10214bbf5ea8SMatthew G. Knepley }
10224bbf5ea8SMatthew G. Knepley 
10234bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchZeroFillMatrix_Private(Mat mat, const PetscInt ncell, const PetscInt ndof, const PetscInt *dof)
10244bbf5ea8SMatthew G. Knepley {
10254bbf5ea8SMatthew G. Knepley   const PetscScalar *values = NULL;
10264bbf5ea8SMatthew G. Knepley   PetscInt           rows, c, i;
10274bbf5ea8SMatthew G. Knepley   PetscErrorCode     ierr;
10284bbf5ea8SMatthew G. Knepley 
10294bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
10304bbf5ea8SMatthew G. Knepley   ierr = PetscCalloc1(ndof*ndof, &values);CHKERRQ(ierr);
10314bbf5ea8SMatthew G. Knepley   for (c = 0; c < ncell; ++c) {
10324bbf5ea8SMatthew G. Knepley     const PetscInt *idx = &dof[ndof*c];
10334bbf5ea8SMatthew G. Knepley     ierr = MatSetValues(mat, ndof, idx, ndof, idx, values, INSERT_VALUES);CHKERRQ(ierr);
10344bbf5ea8SMatthew G. Knepley   }
10354bbf5ea8SMatthew G. Knepley   ierr = MatGetLocalSize(mat, &rows, NULL);CHKERRQ(ierr);
10364bbf5ea8SMatthew G. Knepley   for (i = 0; i < rows; ++i) {
10374bbf5ea8SMatthew G. Knepley     ierr = MatSetValues(mat, 1, &i, 1, &i, values, INSERT_VALUES);CHKERRQ(ierr);
10384bbf5ea8SMatthew G. Knepley   }
10394bbf5ea8SMatthew G. Knepley   ierr = MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10404bbf5ea8SMatthew G. Knepley   ierr = MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10414bbf5ea8SMatthew G. Knepley   ierr = PetscFree(values);CHKERRQ(ierr);
10424bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
10434bbf5ea8SMatthew G. Knepley }
10444bbf5ea8SMatthew G. Knepley 
10454bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchCreateMatrix_Private(PC pc, PetscInt point, Mat *mat)
10464bbf5ea8SMatthew G. Knepley {
10474bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
10484bbf5ea8SMatthew G. Knepley   Vec            x, y;
10494bbf5ea8SMatthew G. Knepley   PetscBool      flg;
10504bbf5ea8SMatthew G. Knepley   PetscInt       csize, rsize;
10514bbf5ea8SMatthew G. Knepley   const char    *prefix = NULL;
10524bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
10534bbf5ea8SMatthew G. Knepley 
10544bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
10554bbf5ea8SMatthew G. Knepley   x = patch->patchX[point];
10564bbf5ea8SMatthew G. Knepley   y = patch->patchY[point];
10574bbf5ea8SMatthew G. Knepley   ierr = VecGetSize(x, &csize);CHKERRQ(ierr);
10584bbf5ea8SMatthew G. Knepley   ierr = VecGetSize(y, &rsize);CHKERRQ(ierr);
10594bbf5ea8SMatthew G. Knepley   ierr = MatCreate(PETSC_COMM_SELF, mat);CHKERRQ(ierr);
10604bbf5ea8SMatthew G. Knepley   ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
10614bbf5ea8SMatthew G. Knepley   ierr = MatSetOptionsPrefix(*mat, prefix);CHKERRQ(ierr);
10624bbf5ea8SMatthew G. Knepley   ierr = MatAppendOptionsPrefix(*mat, "sub_");CHKERRQ(ierr);
10634bbf5ea8SMatthew G. Knepley   if (patch->sub_mat_type) {ierr = MatSetType(*mat, patch->sub_mat_type);CHKERRQ(ierr);}
10644bbf5ea8SMatthew G. Knepley   ierr = MatSetSizes(*mat, rsize, csize, rsize, csize);CHKERRQ(ierr);
10654bbf5ea8SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) *mat, MATDENSE, &flg);CHKERRQ(ierr);
10664bbf5ea8SMatthew G. Knepley   if (!flg) {ierr = PetscObjectTypeCompare((PetscObject)*mat, MATSEQDENSE, &flg);CHKERRQ(ierr);}
10674bbf5ea8SMatthew G. Knepley   /* Sparse patch matrices */
10684bbf5ea8SMatthew G. Knepley   if (!flg) {
10694bbf5ea8SMatthew G. Knepley     PetscBT         bt;
10704bbf5ea8SMatthew G. Knepley     PetscInt       *dnnz      = NULL;
10714bbf5ea8SMatthew G. Knepley     const PetscInt *dofsArray = NULL;
10724bbf5ea8SMatthew G. Knepley     PetscInt        pStart, pEnd, ncell, offset, c, i, j;
10734bbf5ea8SMatthew G. Knepley 
10744bbf5ea8SMatthew G. Knepley     ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
10754bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
10764bbf5ea8SMatthew G. Knepley     point += pStart;
10774bbf5ea8SMatthew G. Knepley     if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
10784bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
10794bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
10804bbf5ea8SMatthew G. Knepley     ierr = PetscCalloc1(rsize, &dnnz);CHKERRQ(ierr);
10814bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventBegin(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
10824bbf5ea8SMatthew G. Knepley     /* XXX: This uses N^2 bits to store the sparsity pattern on a
10834bbf5ea8SMatthew G. Knepley      * patch.  This is probably OK if the patches are not too big,
10844bbf5ea8SMatthew G. Knepley      * but could use quite a bit of memory for planes in 3D.
10854bbf5ea8SMatthew G. Knepley      * Should we switch based on the value of rsize to a
10864bbf5ea8SMatthew G. Knepley      * hash-table (slower, but more memory efficient) approach? */
10874bbf5ea8SMatthew G. Knepley     ierr = PetscBTCreate(rsize*rsize, &bt);CHKERRQ(ierr);
10884bbf5ea8SMatthew G. Knepley     for (c = 0; c < ncell; ++c) {
10894bbf5ea8SMatthew G. Knepley       const PetscInt *idx = dofsArray + (offset + c)*patch->totalDofsPerCell;
10904bbf5ea8SMatthew G. Knepley       for (i = 0; i < patch->totalDofsPerCell; ++i) {
10914bbf5ea8SMatthew G. Knepley         const PetscInt row = idx[i];
10924bbf5ea8SMatthew G. Knepley         for (j = 0; j < patch->totalDofsPerCell; ++j) {
10934bbf5ea8SMatthew G. Knepley           const PetscInt col = idx[j];
10944bbf5ea8SMatthew G. Knepley           const PetscInt key = row*rsize + col;
10954bbf5ea8SMatthew G. Knepley           if (!PetscBTLookupSet(bt, key)) ++dnnz[row];
10964bbf5ea8SMatthew G. Knepley         }
10974bbf5ea8SMatthew G. Knepley       }
10984bbf5ea8SMatthew G. Knepley     }
10994bbf5ea8SMatthew G. Knepley     ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
11004bbf5ea8SMatthew G. Knepley     ierr = MatXAIJSetPreallocation(*mat, 1, dnnz, NULL, NULL, NULL);CHKERRQ(ierr);
11014bbf5ea8SMatthew G. Knepley     ierr = PetscFree(dnnz);CHKERRQ(ierr);
11024bbf5ea8SMatthew G. Knepley     ierr = PCPatchZeroFillMatrix_Private(*mat, ncell, patch->totalDofsPerCell, &dofsArray[offset*patch->totalDofsPerCell]);CHKERRQ(ierr);
11034bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventEnd(PC_Patch_Prealloc, pc, 0, 0, 0);CHKERRQ(ierr);
11044bbf5ea8SMatthew G. Knepley     ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
11054bbf5ea8SMatthew G. Knepley   }
11064bbf5ea8SMatthew G. Knepley   ierr = MatSetUp(*mat);CHKERRQ(ierr);
11074bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
11084bbf5ea8SMatthew G. Knepley }
11094bbf5ea8SMatthew G. Knepley 
11104bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatchComputeOperator_Private(PC pc, Mat mat, Mat multMat, PetscInt point)
11114bbf5ea8SMatthew G. Knepley {
11124bbf5ea8SMatthew G. Knepley   PC_PATCH       *patch = (PC_PATCH *) pc->data;
11134bbf5ea8SMatthew G. Knepley   const PetscInt *dofsArray;
11144bbf5ea8SMatthew G. Knepley   const PetscInt *cellsArray;
11154bbf5ea8SMatthew G. Knepley   PetscInt        ncell, offset, pStart, pEnd;
11164bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
11174bbf5ea8SMatthew G. Knepley 
11184bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
11194bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventBegin(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
11204bbf5ea8SMatthew G. Knepley   if (!patch->usercomputeop) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call PCPatchSetComputeOperator() to set user callback\n");
11214bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
11224bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
11234bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetChart(patch->cellCounts, &pStart, &pEnd);CHKERRQ(ierr);
11244bbf5ea8SMatthew G. Knepley 
11254bbf5ea8SMatthew G. Knepley   point += pStart;
11264bbf5ea8SMatthew G. Knepley   if (point >= pEnd) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Operator point %D not in [%D, %D)\n", point, pStart, pEnd);CHKERRQ(ierr);
11274bbf5ea8SMatthew G. Knepley 
11284bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetDof(patch->cellCounts, point, &ncell);CHKERRQ(ierr);
11294bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetOffset(patch->cellCounts, point, &offset);CHKERRQ(ierr);
11304bbf5ea8SMatthew G. Knepley   if (ncell <= 0) {
11314bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
11324bbf5ea8SMatthew G. Knepley     PetscFunctionReturn(0);
11334bbf5ea8SMatthew G. Knepley   }
11344bbf5ea8SMatthew G. Knepley   PetscStackPush("PCPatch user callback");
11354bbf5ea8SMatthew G. Knepley   ierr = patch->usercomputeop(pc, mat, ncell, cellsArray + offset, ncell*patch->totalDofsPerCell, dofsArray + offset*patch->totalDofsPerCell, patch->usercomputectx);CHKERRQ(ierr);
11364bbf5ea8SMatthew G. Knepley   PetscStackPop;
11374bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(patch->dofs, &dofsArray);CHKERRQ(ierr);
11384bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(patch->cells, &cellsArray);CHKERRQ(ierr);
11394bbf5ea8SMatthew G. Knepley   /* Apply boundary conditions.  Could also do this through the local_to_patch guy. */
11404bbf5ea8SMatthew G. Knepley   if (patch->multiplicative) {
11414bbf5ea8SMatthew G. Knepley     ierr = MatCopy(mat, multMat, SAME_NONZERO_PATTERN);CHKERRQ(ierr);
11424bbf5ea8SMatthew G. Knepley     ierr = MatZeroRowsColumnsIS(multMat, patch->multBcs[point-pStart], (PetscScalar) 1.0, NULL, NULL);CHKERRQ(ierr);
11434bbf5ea8SMatthew G. Knepley   }
11444bbf5ea8SMatthew G. Knepley   ierr = MatZeroRowsColumnsIS(mat, patch->bcs[point-pStart], (PetscScalar) 1.0, NULL, NULL);CHKERRQ(ierr);
11454bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventEnd(PC_Patch_ComputeOp, pc, 0, 0, 0);CHKERRQ(ierr);
11464bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
11474bbf5ea8SMatthew G. Knepley }
11484bbf5ea8SMatthew G. Knepley 
11494bbf5ea8SMatthew G. Knepley static PetscErrorCode PCPatch_ScatterLocal_Private(PC pc, PetscInt p, Vec x, Vec y, InsertMode mode, ScatterMode scat)
11504bbf5ea8SMatthew G. Knepley {
11514bbf5ea8SMatthew G. Knepley   PC_PATCH          *patch     = (PC_PATCH *) pc->data;
11524bbf5ea8SMatthew G. Knepley   const PetscScalar *xArray    = NULL;
11534bbf5ea8SMatthew G. Knepley   PetscScalar       *yArray    = NULL;
11544bbf5ea8SMatthew G. Knepley   const PetscInt    *gtolArray = NULL;
11554bbf5ea8SMatthew G. Knepley   PetscInt           dof, offset, lidx;
11564bbf5ea8SMatthew G. Knepley   PetscErrorCode     ierr;
11574bbf5ea8SMatthew G. Knepley 
11584bbf5ea8SMatthew G. Knepley   PetscFunctionBeginHot;
11594bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventBegin(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
11604bbf5ea8SMatthew G. Knepley   ierr = VecGetArrayRead(x, &xArray);CHKERRQ(ierr);
11614bbf5ea8SMatthew G. Knepley   ierr = VecGetArray(y, &yArray);CHKERRQ(ierr);
11624bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
11634bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetOffset(patch->gtolCounts, p, &offset);CHKERRQ(ierr);
11644bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
11654bbf5ea8SMatthew G. Knepley   if (mode == INSERT_VALUES && scat != SCATTER_FORWARD) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't insert if not scattering forward\n");
11664bbf5ea8SMatthew G. Knepley   if (mode == ADD_VALUES    && scat != SCATTER_REVERSE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Can't add if not scattering reverse\n");
11674bbf5ea8SMatthew G. Knepley   for (lidx = 0; lidx < dof; ++lidx) {
11684bbf5ea8SMatthew G. Knepley     const PetscInt gidx = gtolArray[offset+lidx];
11694bbf5ea8SMatthew G. Knepley 
11704bbf5ea8SMatthew G. Knepley     if (mode == INSERT_VALUES) yArray[lidx]  = xArray[gidx]; /* Forward */
11714bbf5ea8SMatthew G. Knepley     else                       yArray[gidx] += xArray[lidx]; /* Reverse */
11724bbf5ea8SMatthew G. Knepley   }
11734bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(patch->gtol, &gtolArray);CHKERRQ(ierr);
11744bbf5ea8SMatthew G. Knepley   ierr = VecRestoreArrayRead(x, &xArray);CHKERRQ(ierr);
11754bbf5ea8SMatthew G. Knepley   ierr = VecRestoreArray(y, &yArray);CHKERRQ(ierr);
11764bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventEnd(PC_Patch_Scatter, pc, 0, 0, 0);CHKERRQ(ierr);
11774bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
11784bbf5ea8SMatthew G. Knepley }
11794bbf5ea8SMatthew G. Knepley 
11804bbf5ea8SMatthew G. Knepley static PetscErrorCode PCSetUp_PATCH(PC pc)
11814bbf5ea8SMatthew G. Knepley {
11824bbf5ea8SMatthew G. Knepley   PC_PATCH       *patch   = (PC_PATCH *) pc->data;
11834bbf5ea8SMatthew G. Knepley   PetscScalar    *patchX  = NULL;
11844bbf5ea8SMatthew G. Knepley   const PetscInt *bcNodes = NULL;
11854bbf5ea8SMatthew G. Knepley   PetscInt        numBcs, i, j;
11864bbf5ea8SMatthew G. Knepley   const char     *prefix;
11874bbf5ea8SMatthew G. Knepley   PetscErrorCode  ierr;
11884bbf5ea8SMatthew G. Knepley 
11894bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
11904bbf5ea8SMatthew G. Knepley   if (!pc->setupcalled) {
11914bbf5ea8SMatthew G. Knepley     PetscInt pStart, pEnd, p;
11924bbf5ea8SMatthew G. Knepley     PetscInt localSize;
11934bbf5ea8SMatthew G. Knepley 
11944bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventBegin(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
11954bbf5ea8SMatthew G. Knepley 
11964bbf5ea8SMatthew G. Knepley     localSize = patch->subspaceOffsets[patch->nsubspaces];
11974bbf5ea8SMatthew G. Knepley     ierr = VecCreateSeq(PETSC_COMM_SELF, localSize, &patch->localX);CHKERRQ(ierr);
11984bbf5ea8SMatthew G. Knepley     ierr = VecSetUp(patch->localX);CHKERRQ(ierr);
11994bbf5ea8SMatthew G. Knepley     ierr = VecDuplicate(patch->localX, &patch->localY);CHKERRQ(ierr);
12004bbf5ea8SMatthew G. Knepley     ierr = PCPatchCreateCellPatches(pc);CHKERRQ(ierr);
12014bbf5ea8SMatthew G. Knepley     ierr = PCPatchCreateCellPatchDiscretisationInfo(pc);CHKERRQ(ierr);
12024bbf5ea8SMatthew G. Knepley     ierr = PCPatchCreateCellPatchBCs(pc);CHKERRQ(ierr);
12034bbf5ea8SMatthew G. Knepley 
12044bbf5ea8SMatthew G. Knepley     /* OK, now build the work vectors */
12054bbf5ea8SMatthew G. Knepley     ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, &pEnd);CHKERRQ(ierr);
12064bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(patch->npatch, &patch->patchX);CHKERRQ(ierr);
12074bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(patch->npatch, &patch->patchY);CHKERRQ(ierr);
12084bbf5ea8SMatthew G. Knepley     if (patch->partition_of_unity && patch->multiplicative) {ierr = PetscMalloc1(patch->npatch, &patch->patch_dof_weights);CHKERRQ(ierr);}
12094bbf5ea8SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
12104bbf5ea8SMatthew G. Knepley       PetscInt dof;
12114bbf5ea8SMatthew G. Knepley 
12124bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetDof(patch->gtolCounts, p, &dof);CHKERRQ(ierr);
12134bbf5ea8SMatthew G. Knepley       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchX[p-pStart]);CHKERRQ(ierr);
12144bbf5ea8SMatthew G. Knepley       ierr = VecSetUp(patch->patchX[p-pStart]);CHKERRQ(ierr);
12154bbf5ea8SMatthew G. Knepley       ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patchY[p-pStart]);CHKERRQ(ierr);
12164bbf5ea8SMatthew G. Knepley       ierr = VecSetUp(patch->patchY[p-pStart]);CHKERRQ(ierr);
12174bbf5ea8SMatthew G. Knepley       if (patch->partition_of_unity && patch->multiplicative) {
12184bbf5ea8SMatthew G. Knepley         ierr = VecCreateSeq(PETSC_COMM_SELF, dof, &patch->patch_dof_weights[p-pStart]);CHKERRQ(ierr);
12194bbf5ea8SMatthew G. Knepley         ierr = VecSetUp(patch->patch_dof_weights[p-pStart]);CHKERRQ(ierr);
12204bbf5ea8SMatthew G. Knepley       }
12214bbf5ea8SMatthew G. Knepley     }
12224bbf5ea8SMatthew G. Knepley     ierr = PetscMalloc1(patch->npatch, &patch->ksp);CHKERRQ(ierr);
12234bbf5ea8SMatthew G. Knepley     ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
12244bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
12254bbf5ea8SMatthew G. Knepley       ierr = KSPCreate(PETSC_COMM_SELF, &patch->ksp[i]);CHKERRQ(ierr);
12264bbf5ea8SMatthew G. Knepley       ierr = KSPSetOptionsPrefix(patch->ksp[i], prefix);CHKERRQ(ierr);
12274bbf5ea8SMatthew G. Knepley       ierr = KSPAppendOptionsPrefix(patch->ksp[i], "sub_");CHKERRQ(ierr);
12284bbf5ea8SMatthew G. Knepley     }
12294bbf5ea8SMatthew G. Knepley     if (patch->save_operators) {
12304bbf5ea8SMatthew G. Knepley       ierr = PetscMalloc1(patch->npatch, &patch->mat);CHKERRQ(ierr);
12314bbf5ea8SMatthew G. Knepley       if (patch->multiplicative) {ierr = PetscMalloc1(patch->npatch, &patch->multMat);CHKERRQ(ierr);}
12324bbf5ea8SMatthew G. Knepley       for (i = 0; i < patch->npatch; ++i) {
12334bbf5ea8SMatthew G. Knepley         ierr = PCPatchCreateMatrix_Private(pc, i, &patch->mat[i]);CHKERRQ(ierr);
12344bbf5ea8SMatthew G. Knepley         if (patch->multiplicative) {ierr = MatDuplicate(patch->mat[i], MAT_SHARE_NONZERO_PATTERN, &patch->multMat[i]);CHKERRQ(ierr);}
12354bbf5ea8SMatthew G. Knepley       }
12364bbf5ea8SMatthew G. Knepley     }
12374bbf5ea8SMatthew G. Knepley     ierr = PetscLogEventEnd(PC_Patch_CreatePatches, pc, 0, 0, 0);CHKERRQ(ierr);
12384bbf5ea8SMatthew G. Knepley 
12394bbf5ea8SMatthew G. Knepley     /* If desired, calculate weights for dof multiplicity */
12404bbf5ea8SMatthew G. Knepley     if (patch->partition_of_unity) {
12414bbf5ea8SMatthew G. Knepley       ierr = VecDuplicate(patch->localX, &patch->dof_weights);CHKERRQ(ierr);
12424bbf5ea8SMatthew G. Knepley       for (i = 0; i < patch->npatch; ++i) {
12434bbf5ea8SMatthew G. Knepley         PetscInt dof;
12444bbf5ea8SMatthew G. Knepley 
12454bbf5ea8SMatthew G. Knepley         ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &dof);CHKERRQ(ierr);
12464bbf5ea8SMatthew G. Knepley         if (dof <= 0) continue;
12474bbf5ea8SMatthew G. Knepley         ierr = VecSet(patch->patchX[i], 1.0);CHKERRQ(ierr);
12484bbf5ea8SMatthew G. Knepley         /* TODO: Do we need different scatters for X and Y? */
12494bbf5ea8SMatthew G. Knepley         ierr = VecGetArray(patch->patchX[i], &patchX);CHKERRQ(ierr);
12504bbf5ea8SMatthew G. Knepley         /* Apply bcs to patchX (zero entries) */
12514bbf5ea8SMatthew G. Knepley         ierr = ISGetLocalSize(patch->bcs[i], &numBcs);CHKERRQ(ierr);
12524bbf5ea8SMatthew G. Knepley         ierr = ISGetIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr);
12534bbf5ea8SMatthew G. Knepley         for (j = 0; j < numBcs; ++j) patchX[bcNodes[j]] = 0;
12544bbf5ea8SMatthew G. Knepley         ierr = ISRestoreIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr);
12554bbf5ea8SMatthew G. Knepley         ierr = VecRestoreArray(patch->patchX[i], &patchX);CHKERRQ(ierr);
12564bbf5ea8SMatthew G. Knepley 
12574bbf5ea8SMatthew G. Knepley         ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchX[i], patch->dof_weights, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
12584bbf5ea8SMatthew G. Knepley       }
12594bbf5ea8SMatthew G. Knepley       ierr = VecReciprocal(patch->dof_weights);CHKERRQ(ierr);
12604bbf5ea8SMatthew G. Knepley       if (patch->partition_of_unity && patch->multiplicative) {
12614bbf5ea8SMatthew G. Knepley         for (i = 0; i < patch->npatch; ++i) {
12624bbf5ea8SMatthew G. Knepley           ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->dof_weights, patch->patch_dof_weights[i], INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
12634bbf5ea8SMatthew G. Knepley         }
12644bbf5ea8SMatthew G. Knepley       }
12654bbf5ea8SMatthew G. Knepley     }
12664bbf5ea8SMatthew G. Knepley   }
12674bbf5ea8SMatthew G. Knepley   if (patch->save_operators) {
12684bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
12694bbf5ea8SMatthew G. Knepley       ierr = MatZeroEntries(patch->mat[i]);CHKERRQ(ierr);
12704bbf5ea8SMatthew G. Knepley       if (patch->multiplicative) {ierr = PCPatchComputeOperator_Private(pc, patch->mat[i], patch->multMat[i], i);CHKERRQ(ierr);}
12714bbf5ea8SMatthew G. Knepley       else                       {ierr = PCPatchComputeOperator_Private(pc, patch->mat[i], NULL, i);CHKERRQ(ierr);}
12724bbf5ea8SMatthew G. Knepley       ierr = KSPSetOperators(patch->ksp[i], patch->mat[i], patch->mat[i]);CHKERRQ(ierr);
12734bbf5ea8SMatthew G. Knepley     }
12744bbf5ea8SMatthew G. Knepley   }
12754bbf5ea8SMatthew G. Knepley   if (!pc->setupcalled) {
12764bbf5ea8SMatthew G. Knepley     /* TODO: Only call if SetFromOptions was called on this PC */
12774bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {ierr = KSPSetFromOptions(patch->ksp[i]);CHKERRQ(ierr);}
12784bbf5ea8SMatthew G. Knepley   }
12794bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
12804bbf5ea8SMatthew G. Knepley }
12814bbf5ea8SMatthew G. Knepley 
12824bbf5ea8SMatthew G. Knepley static PetscErrorCode PCApply_PATCH(PC pc, Vec x, Vec y)
12834bbf5ea8SMatthew G. Knepley {
12844bbf5ea8SMatthew G. Knepley   PC_PATCH          *patch    = (PC_PATCH *) pc->data;
12854bbf5ea8SMatthew G. Knepley   const PetscScalar *globalX  = NULL;
12864bbf5ea8SMatthew G. Knepley   PetscScalar       *localX   = NULL;
12874bbf5ea8SMatthew G. Knepley   PetscScalar       *globalY  = NULL;
12884bbf5ea8SMatthew G. Knepley   PetscScalar       *patchX   = NULL;
12894bbf5ea8SMatthew G. Knepley   const PetscInt    *bcNodes  = NULL;
12904bbf5ea8SMatthew G. Knepley   PetscInt           nsweep   = patch->symmetrise_sweep ? 2 : 1;
12914bbf5ea8SMatthew G. Knepley   PetscInt           start[2] = {0, 0};
12924bbf5ea8SMatthew G. Knepley   PetscInt           end[2]   = {-1, -1};
12934bbf5ea8SMatthew G. Knepley   const PetscInt     inc[2]   = {1, -1};
12944bbf5ea8SMatthew G. Knepley   const PetscScalar *localY;
12954bbf5ea8SMatthew G. Knepley   const PetscInt    *iterationSet;
12964bbf5ea8SMatthew G. Knepley   PetscInt           pStart, numBcs, n, sweep, bc, j;
12974bbf5ea8SMatthew G. Knepley   PetscErrorCode     ierr;
12984bbf5ea8SMatthew G. Knepley 
12994bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
13004bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventBegin(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
13014bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsPushGetViewerOff(PETSC_TRUE);CHKERRQ(ierr);
13024bbf5ea8SMatthew G. Knepley   end[0]   = patch->npatch;
13034bbf5ea8SMatthew G. Knepley   start[1] = patch->npatch-1;
13044bbf5ea8SMatthew G. Knepley   if (patch->user_patches) {
13054bbf5ea8SMatthew G. Knepley     ierr = ISGetLocalSize(patch->iterationSet, &end[0]);CHKERRQ(ierr);
13064bbf5ea8SMatthew G. Knepley     start[1] = end[0] - 1;
13074bbf5ea8SMatthew G. Knepley     ierr = ISGetIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);
13084bbf5ea8SMatthew G. Knepley   }
13094bbf5ea8SMatthew G. Knepley   /* Scatter from global space into overlapped local spaces */
13104bbf5ea8SMatthew G. Knepley   ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr);
13114bbf5ea8SMatthew G. Knepley   ierr = VecGetArray(patch->localX, &localX);CHKERRQ(ierr);
13124bbf5ea8SMatthew G. Knepley   ierr = PetscSFBcastBegin(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr);
13134bbf5ea8SMatthew G. Knepley   ierr = PetscSFBcastEnd(patch->defaultSF, MPIU_SCALAR, globalX, localX);CHKERRQ(ierr);
13144bbf5ea8SMatthew G. Knepley   ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr);
13154bbf5ea8SMatthew G. Knepley   ierr = VecRestoreArray(patch->localX, &localX);CHKERRQ(ierr);
13164bbf5ea8SMatthew G. Knepley 
13174bbf5ea8SMatthew G. Knepley   ierr = VecSet(patch->localY, 0.0);CHKERRQ(ierr);
13184bbf5ea8SMatthew G. Knepley   ierr = PetscSectionGetChart(patch->gtolCounts, &pStart, NULL);CHKERRQ(ierr);
13194bbf5ea8SMatthew G. Knepley   for (sweep = 0; sweep < nsweep; sweep++) {
13204bbf5ea8SMatthew G. Knepley     for (j = start[sweep]; j*inc[sweep] < end[sweep]*inc[sweep]; j += inc[sweep]) {
13214bbf5ea8SMatthew G. Knepley       Mat      multMat = NULL;
13224bbf5ea8SMatthew G. Knepley       PetscInt i       = patch->user_patches ? iterationSet[j] : j;
13234bbf5ea8SMatthew G. Knepley       PetscInt start, len;
13244bbf5ea8SMatthew G. Knepley 
13254bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetDof(patch->gtolCounts, i+pStart, &len);CHKERRQ(ierr);
13264bbf5ea8SMatthew G. Knepley       ierr = PetscSectionGetOffset(patch->gtolCounts, i+pStart, &start);CHKERRQ(ierr);
13274bbf5ea8SMatthew G. Knepley       /* TODO: Squash out these guys in the setup as well. */
13284bbf5ea8SMatthew G. Knepley       if (len <= 0) continue;
13294bbf5ea8SMatthew G. Knepley       /* TODO: Do we need different scatters for X and Y? */
13304bbf5ea8SMatthew G. Knepley       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->localX, patch->patchX[i], INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
13314bbf5ea8SMatthew G. Knepley       /* Apply bcs to patchX (zero entries) */
13324bbf5ea8SMatthew G. Knepley       ierr = VecGetArray(patch->patchX[i], &patchX);CHKERRQ(ierr);
13334bbf5ea8SMatthew G. Knepley       ierr = ISGetLocalSize(patch->bcs[i], &numBcs);CHKERRQ(ierr);
13344bbf5ea8SMatthew G. Knepley       ierr = ISGetIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr);
13354bbf5ea8SMatthew G. Knepley       for (bc = 0; bc < numBcs; ++bc) patchX[bcNodes[bc]] = 0;
13364bbf5ea8SMatthew G. Knepley       ierr = ISRestoreIndices(patch->bcs[i], &bcNodes);CHKERRQ(ierr);
13374bbf5ea8SMatthew G. Knepley       ierr = VecRestoreArray(patch->patchX[i], &patchX);CHKERRQ(ierr);
13384bbf5ea8SMatthew G. Knepley       if (!patch->save_operators) {
13394bbf5ea8SMatthew G. Knepley         Mat mat;
13404bbf5ea8SMatthew G. Knepley 
13414bbf5ea8SMatthew G. Knepley         ierr = PCPatchCreateMatrix_Private(pc, i, &mat);CHKERRQ(ierr);
13424bbf5ea8SMatthew G. Knepley         if (patch->multiplicative) {ierr = MatDuplicate(mat, MAT_SHARE_NONZERO_PATTERN, &multMat);CHKERRQ(ierr);}
13434bbf5ea8SMatthew G. Knepley         /* Populate operator here. */
13444bbf5ea8SMatthew G. Knepley         ierr = PCPatchComputeOperator_Private(pc, mat, multMat, i);CHKERRQ(ierr);
13454bbf5ea8SMatthew G. Knepley         ierr = KSPSetOperators(patch->ksp[i], mat, mat);
13464bbf5ea8SMatthew G. Knepley         /* Drop reference so the KSPSetOperators below will blow it away. */
13474bbf5ea8SMatthew G. Knepley         ierr = MatDestroy(&mat);CHKERRQ(ierr);
13484bbf5ea8SMatthew G. Knepley       } else if (patch->multiplicative) {
13494bbf5ea8SMatthew G. Knepley         multMat = patch->multMat[i];
13504bbf5ea8SMatthew G. Knepley       }
13514bbf5ea8SMatthew G. Knepley 
13524bbf5ea8SMatthew G. Knepley       ierr = PetscLogEventBegin(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
13534bbf5ea8SMatthew G. Knepley       ierr = KSPSolve(patch->ksp[i], patch->patchX[i], patch->patchY[i]);CHKERRQ(ierr);
13544bbf5ea8SMatthew G. Knepley       ierr = PetscLogEventEnd(PC_Patch_Solve, pc, 0, 0, 0);CHKERRQ(ierr);
13554bbf5ea8SMatthew G. Knepley 
13564bbf5ea8SMatthew G. Knepley       if (!patch->save_operators) {
13574bbf5ea8SMatthew G. Knepley         PC pc;
13584bbf5ea8SMatthew G. Knepley         ierr = KSPSetOperators(patch->ksp[i], NULL, NULL);CHKERRQ(ierr);
13594bbf5ea8SMatthew G. Knepley         ierr = KSPGetPC(patch->ksp[i], &pc);CHKERRQ(ierr);
13604bbf5ea8SMatthew G. Knepley         /* Destroy PC context too, otherwise the factored matrix hangs around. */
13614bbf5ea8SMatthew G. Knepley         ierr = PCReset(pc);CHKERRQ(ierr);
13624bbf5ea8SMatthew G. Knepley       }
13634bbf5ea8SMatthew G. Knepley 
13644bbf5ea8SMatthew G. Knepley       if (patch->partition_of_unity && patch->multiplicative) {
13654bbf5ea8SMatthew G. Knepley         ierr = VecPointwiseMult(patch->patchY[i], patch->patchY[i], patch->patch_dof_weights[i]);CHKERRQ(ierr);
13664bbf5ea8SMatthew G. Knepley       }
13674bbf5ea8SMatthew G. Knepley       ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchY[i], patch->localY, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
13684bbf5ea8SMatthew G. Knepley 
13694bbf5ea8SMatthew G. Knepley       /* Get the matrix on the patch but with only global bcs applied.
13704bbf5ea8SMatthew G. Knepley        * This matrix is then multiplied with the result from the previous solve
13714bbf5ea8SMatthew G. Knepley        * to obtain by how much the residual changes. */
13724bbf5ea8SMatthew G. Knepley       if (patch->multiplicative) {
13734bbf5ea8SMatthew G. Knepley         ierr = MatMult(multMat, patch->patchY[i], patch->patchX[i]);CHKERRQ(ierr);
13744bbf5ea8SMatthew G. Knepley         ierr = VecScale(patch->patchX[i], -1.0);CHKERRQ(ierr);
13754bbf5ea8SMatthew G. Knepley         ierr = PCPatch_ScatterLocal_Private(pc, i+pStart, patch->patchX[i], patch->localX, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
13764bbf5ea8SMatthew G. Knepley         if (!patch->save_operators) {ierr = MatDestroy(&multMat);CHKERRQ(ierr);}
13774bbf5ea8SMatthew G. Knepley       }
13784bbf5ea8SMatthew G. Knepley     }
13794bbf5ea8SMatthew G. Knepley   }
13804bbf5ea8SMatthew G. Knepley   if (patch->user_patches) {ierr = ISRestoreIndices(patch->iterationSet, &iterationSet);CHKERRQ(ierr);}
13814bbf5ea8SMatthew G. Knepley   /* XXX: should we do this on the global vector? */
13824bbf5ea8SMatthew G. Knepley   if (patch->partition_of_unity && !patch->multiplicative) {
13834bbf5ea8SMatthew G. Knepley     ierr = VecPointwiseMult(patch->localY, patch->localY, patch->dof_weights);CHKERRQ(ierr);
13844bbf5ea8SMatthew G. Knepley   }
13854bbf5ea8SMatthew G. Knepley   /* Now patch->localY contains the solution of the patch solves, so we need to combine them all. */
13864bbf5ea8SMatthew G. Knepley   ierr = VecSet(y, 0.0);CHKERRQ(ierr);
13874bbf5ea8SMatthew G. Knepley   ierr = VecGetArray(y, &globalY);CHKERRQ(ierr);
13884bbf5ea8SMatthew G. Knepley   ierr = VecGetArrayRead(patch->localY, &localY);CHKERRQ(ierr);
13894bbf5ea8SMatthew G. Knepley   ierr = PetscSFReduceBegin(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr);
13904bbf5ea8SMatthew G. Knepley   ierr = PetscSFReduceEnd(patch->defaultSF, MPIU_SCALAR, localY, globalY, MPI_SUM);CHKERRQ(ierr);
13914bbf5ea8SMatthew G. Knepley   ierr = VecRestoreArrayRead(patch->localY, &localY);CHKERRQ(ierr);
13924bbf5ea8SMatthew G. Knepley 
13934bbf5ea8SMatthew G. Knepley   /* Now we need to send the global BC values through */
13944bbf5ea8SMatthew G. Knepley   ierr = VecGetArrayRead(x, &globalX);CHKERRQ(ierr);
13954bbf5ea8SMatthew G. Knepley   ierr = ISGetSize(patch->globalBcNodes, &numBcs);CHKERRQ(ierr);
13964bbf5ea8SMatthew G. Knepley   ierr = ISGetIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
13974bbf5ea8SMatthew G. Knepley   ierr = VecGetLocalSize(x, &n);CHKERRQ(ierr);
13984bbf5ea8SMatthew G. Knepley   for (bc = 0; bc < numBcs; ++bc) {
13994bbf5ea8SMatthew G. Knepley     const PetscInt idx = bcNodes[bc];
14004bbf5ea8SMatthew G. Knepley     if (idx < n) globalY[idx] = globalX[idx];
14014bbf5ea8SMatthew G. Knepley   }
14024bbf5ea8SMatthew G. Knepley 
14034bbf5ea8SMatthew G. Knepley   ierr = ISRestoreIndices(patch->globalBcNodes, &bcNodes);CHKERRQ(ierr);
14044bbf5ea8SMatthew G. Knepley   ierr = VecRestoreArrayRead(x, &globalX);CHKERRQ(ierr);
14054bbf5ea8SMatthew G. Knepley   ierr = VecRestoreArray(y, &globalY);CHKERRQ(ierr);
14064bbf5ea8SMatthew G. Knepley 
14074bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsPopGetViewerOff();CHKERRQ(ierr);
14084bbf5ea8SMatthew G. Knepley   ierr = PetscLogEventEnd(PC_Patch_Apply, pc, 0, 0, 0);CHKERRQ(ierr);
14094bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
14104bbf5ea8SMatthew G. Knepley }
14114bbf5ea8SMatthew G. Knepley 
14124bbf5ea8SMatthew G. Knepley static PetscErrorCode PCReset_PATCH(PC pc)
14134bbf5ea8SMatthew G. Knepley {
14144bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
14154bbf5ea8SMatthew G. Knepley   PetscInt       i;
14164bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
14174bbf5ea8SMatthew G. Knepley 
14184bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
14194bbf5ea8SMatthew G. Knepley   /* TODO: Get rid of all these ifs */
14204bbf5ea8SMatthew G. Knepley   ierr = PetscSFDestroy(&patch->defaultSF);CHKERRQ(ierr);
14214bbf5ea8SMatthew G. Knepley   ierr = PetscSectionDestroy(&patch->cellCounts);CHKERRQ(ierr);
14224bbf5ea8SMatthew G. Knepley   ierr = PetscSectionDestroy(&patch->cellNumbering);CHKERRQ(ierr);
14234bbf5ea8SMatthew G. Knepley   ierr = PetscSectionDestroy(&patch->gtolCounts);CHKERRQ(ierr);
14244bbf5ea8SMatthew G. Knepley   ierr = PetscSectionDestroy(&patch->bcCounts);CHKERRQ(ierr);
14254bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->gtol);CHKERRQ(ierr);
14264bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->cells);CHKERRQ(ierr);
14274bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->dofs);CHKERRQ(ierr);
14284bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->ghostBcNodes);CHKERRQ(ierr);
14294bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->globalBcNodes);CHKERRQ(ierr);
14304bbf5ea8SMatthew G. Knepley 
14314bbf5ea8SMatthew G. Knepley   if (patch->dofSection) {
14324bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->nsubspaces; i++) {
14334bbf5ea8SMatthew G. Knepley       ierr = PetscSectionDestroy(&patch->dofSection[i]);CHKERRQ(ierr);
14344bbf5ea8SMatthew G. Knepley     }
14354bbf5ea8SMatthew G. Knepley   }
14364bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->dofSection);CHKERRQ(ierr);
14374bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->bs);CHKERRQ(ierr);
14384bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->nodesPerCell);CHKERRQ(ierr);
14394bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->cellNodeMap);CHKERRQ(ierr);
14404bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->subspaceOffsets);CHKERRQ(ierr);
14414bbf5ea8SMatthew G. Knepley 
14424bbf5ea8SMatthew G. Knepley   if (patch->bcs) {
14434bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
14444bbf5ea8SMatthew G. Knepley       ierr = ISDestroy(&patch->bcs[i]);CHKERRQ(ierr);
14454bbf5ea8SMatthew G. Knepley     }
14464bbf5ea8SMatthew G. Knepley     ierr = PetscFree(patch->bcs);CHKERRQ(ierr);
14474bbf5ea8SMatthew G. Knepley   }
14484bbf5ea8SMatthew G. Knepley   if (patch->multBcs) {
14494bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
14504bbf5ea8SMatthew G. Knepley       ierr = ISDestroy(&patch->multBcs[i]);CHKERRQ(ierr);
14514bbf5ea8SMatthew G. Knepley     }
14524bbf5ea8SMatthew G. Knepley     ierr = PetscFree(patch->multBcs);CHKERRQ(ierr);
14534bbf5ea8SMatthew G. Knepley   }
14544bbf5ea8SMatthew G. Knepley   if (patch->ksp) {
14554bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
14564bbf5ea8SMatthew G. Knepley       ierr = KSPReset(patch->ksp[i]);CHKERRQ(ierr);
14574bbf5ea8SMatthew G. Knepley     }
14584bbf5ea8SMatthew G. Knepley   }
14594bbf5ea8SMatthew G. Knepley 
14604bbf5ea8SMatthew G. Knepley   ierr = VecDestroy(&patch->localX);CHKERRQ(ierr);
14614bbf5ea8SMatthew G. Knepley   ierr = VecDestroy(&patch->localY);CHKERRQ(ierr);
14624bbf5ea8SMatthew G. Knepley   if (patch->patchX) {
14634bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
14644bbf5ea8SMatthew G. Knepley       ierr = VecDestroy(&patch->patchX[i]);CHKERRQ(ierr);
14654bbf5ea8SMatthew G. Knepley     }
14664bbf5ea8SMatthew G. Knepley     ierr = PetscFree(patch->patchX);CHKERRQ(ierr);
14674bbf5ea8SMatthew G. Knepley   }
14684bbf5ea8SMatthew G. Knepley   if (patch->patchY) {
14694bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
14704bbf5ea8SMatthew G. Knepley       ierr = VecDestroy(&patch->patchY[i]);CHKERRQ(ierr);
14714bbf5ea8SMatthew G. Knepley     }
14724bbf5ea8SMatthew G. Knepley     ierr = PetscFree(patch->patchY);CHKERRQ(ierr);
14734bbf5ea8SMatthew G. Knepley   }
14744bbf5ea8SMatthew G. Knepley   if (patch->partition_of_unity) {
14754bbf5ea8SMatthew G. Knepley     ierr = VecDestroy(&patch->dof_weights);CHKERRQ(ierr);
14764bbf5ea8SMatthew G. Knepley   }
14774bbf5ea8SMatthew G. Knepley   if (patch->patch_dof_weights) {
14784bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
14794bbf5ea8SMatthew G. Knepley       ierr = VecDestroy(&patch->patch_dof_weights[i]);CHKERRQ(ierr);
14804bbf5ea8SMatthew G. Knepley     }
14814bbf5ea8SMatthew G. Knepley     ierr = PetscFree(patch->patch_dof_weights);CHKERRQ(ierr);
14824bbf5ea8SMatthew G. Knepley   }
14834bbf5ea8SMatthew G. Knepley   if (patch->mat) {
14844bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {
14854bbf5ea8SMatthew G. Knepley       ierr = MatDestroy(&patch->mat[i]);CHKERRQ(ierr);
14864bbf5ea8SMatthew G. Knepley       if (patch->multiplicative) {
14874bbf5ea8SMatthew G. Knepley         ierr = MatDestroy(&patch->multMat[i]);CHKERRQ(ierr);
14884bbf5ea8SMatthew G. Knepley       }
14894bbf5ea8SMatthew G. Knepley     }
14904bbf5ea8SMatthew G. Knepley     ierr = PetscFree(patch->mat);CHKERRQ(ierr);
14914bbf5ea8SMatthew G. Knepley     if (patch->multiplicative) {
14924bbf5ea8SMatthew G. Knepley       ierr = PetscFree(patch->multMat);CHKERRQ(ierr);
14934bbf5ea8SMatthew G. Knepley     }
14944bbf5ea8SMatthew G. Knepley   }
14954bbf5ea8SMatthew G. Knepley   ierr = PetscFree(patch->sub_mat_type);CHKERRQ(ierr);
14964bbf5ea8SMatthew G. Knepley 
14974bbf5ea8SMatthew G. Knepley   patch->bs = 0;
14984bbf5ea8SMatthew G. Knepley   patch->cellNodeMap = NULL;
14994bbf5ea8SMatthew G. Knepley 
15004bbf5ea8SMatthew G. Knepley   if (patch->user_patches) {
15014bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->nuserIS; ++i) {
15024bbf5ea8SMatthew G. Knepley       ierr = ISDestroy(&patch->userIS[i]);CHKERRQ(ierr);
15034bbf5ea8SMatthew G. Knepley     }
15044bbf5ea8SMatthew G. Knepley     PetscFree(patch->userIS);
15054bbf5ea8SMatthew G. Knepley     patch->nuserIS = 0;
15064bbf5ea8SMatthew G. Knepley   }
15074bbf5ea8SMatthew G. Knepley   ierr = ISDestroy(&patch->iterationSet);CHKERRQ(ierr);
15084bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
15094bbf5ea8SMatthew G. Knepley }
15104bbf5ea8SMatthew G. Knepley 
15114bbf5ea8SMatthew G. Knepley static PetscErrorCode PCDestroy_PATCH(PC pc)
15124bbf5ea8SMatthew G. Knepley {
15134bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
15144bbf5ea8SMatthew G. Knepley   PetscInt       i;
15154bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
15164bbf5ea8SMatthew G. Knepley 
15174bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
15184bbf5ea8SMatthew G. Knepley   ierr = PCReset_PATCH(pc);CHKERRQ(ierr);
15194bbf5ea8SMatthew G. Knepley   if (patch->ksp) {
15204bbf5ea8SMatthew G. Knepley     for (i = 0; i < patch->npatch; ++i) {ierr = KSPDestroy(&patch->ksp[i]);CHKERRQ(ierr);}
15214bbf5ea8SMatthew G. Knepley     ierr = PetscFree(patch->ksp);CHKERRQ(ierr);
15224bbf5ea8SMatthew G. Knepley   }
15234bbf5ea8SMatthew G. Knepley   ierr = PetscFree(pc->data);CHKERRQ(ierr);
15244bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
15254bbf5ea8SMatthew G. Knepley }
15264bbf5ea8SMatthew G. Knepley 
15274bbf5ea8SMatthew G. Knepley static PetscErrorCode PCSetFromOptions_PATCH(PetscOptionItems *PetscOptionsObject, PC pc)
15284bbf5ea8SMatthew G. Knepley {
15294bbf5ea8SMatthew G. Knepley   PC_PATCH            *patch = (PC_PATCH *) pc->data;
15304bbf5ea8SMatthew G. Knepley   PCPatchConstructType patchConstructionType = PC_PATCH_STAR;
15314bbf5ea8SMatthew G. Knepley   char                 sub_mat_type[256];
15324bbf5ea8SMatthew G. Knepley   PetscBool            flg, dimflg, codimflg;
15334bbf5ea8SMatthew G. Knepley   PetscErrorCode       ierr;
15344bbf5ea8SMatthew G. Knepley 
15354bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
15364bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsHead(PetscOptionsObject, "Vertex-patch Additive Schwarz options");CHKERRQ(ierr);
15374bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_patch_save_operators",  "Store all patch operators for lifetime of PC?", "PCPatchSetSaveOperators", patch->save_operators, &patch->save_operators, &flg);CHKERRQ(ierr);
15384bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_patch_partition_of_unity", "Weight contributions by dof multiplicity?", "PCPatchSetPartitionOfUnity", patch->partition_of_unity, &patch->partition_of_unity, &flg);CHKERRQ(ierr);
15394bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_patch_multiplicative", "Gauss-Seidel instead of Jacobi?", "PCPatchSetMultiplicative", patch->multiplicative, &patch->multiplicative, &flg);CHKERRQ(ierr);
15404bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_patch_construct_dim", "What dimension of mesh point to construct patches by? (0 = vertices)", "PCSetFromOptions_PATCH", patch->dim, &patch->dim, &dimflg);CHKERRQ(ierr);
15414bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_patch_construct_codim", "What co-dimension of mesh point to construct patches by? (0 = cells)", "PCSetFromOptions_PATCH", patch->codim, &patch->codim, &codimflg);CHKERRQ(ierr);
15424bbf5ea8SMatthew G. Knepley   if (dimflg && codimflg) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Can only set one of dimension or co-dimension");CHKERRQ(ierr);
15434bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_patch_construct_type", "How should the patches be constructed?", "PCPatchSetConstructType", PCPatchConstructTypes, (PetscEnum) patchConstructionType, (PetscEnum *) &patchConstructionType, &flg);CHKERRQ(ierr);
15444bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_patch_vanka_dim", "Topological dimension of entities for Vanka to ignore", "PCSetFromOptions_PATCH", patch->vankadim, &patch->vankadim, &flg);
15454bbf5ea8SMatthew G. Knepley   if (flg) {ierr = PCPatchSetConstructType(pc, patchConstructionType, NULL, NULL);CHKERRQ(ierr);}
15464bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsFList("-pc_patch_sub_mat_type", "Matrix type for patch solves", "PCPatchSetSubMatType", MatList, NULL, sub_mat_type, 256, &flg);CHKERRQ(ierr);
15474bbf5ea8SMatthew G. Knepley   if (flg) {ierr = PCPatchSetSubMatType(pc, sub_mat_type);CHKERRQ(ierr);}
15484bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_patch_print_patches", "Print out information during patch construction?", "PCSetFromOptions_PATCH", patch->print_patches, &patch->print_patches, &flg);CHKERRQ(ierr);
15494bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_patch_symmetrise_sweep", "Go start->end, end->start?", "PCSetFromOptions_PATCH", patch->symmetrise_sweep, &patch->symmetrise_sweep, &flg);CHKERRQ(ierr);
15504bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_patch_exclude_subspace", "What subspace (if any) to exclude in construction?", "PCSetFromOptions_PATCH", patch->exclude_subspace, &patch->exclude_subspace, &flg);
15514bbf5ea8SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
15524bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
15534bbf5ea8SMatthew G. Knepley }
15544bbf5ea8SMatthew G. Knepley 
15554bbf5ea8SMatthew G. Knepley static PetscErrorCode PCSetUpOnBlocks_PATCH(PC pc)
15564bbf5ea8SMatthew G. Knepley {
15574bbf5ea8SMatthew G. Knepley   PC_PATCH          *patch = (PC_PATCH*) pc->data;
15584bbf5ea8SMatthew G. Knepley   KSPConvergedReason reason;
15594bbf5ea8SMatthew G. Knepley   PetscInt           i;
15604bbf5ea8SMatthew G. Knepley   PetscErrorCode     ierr;
15614bbf5ea8SMatthew G. Knepley 
15624bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
15634bbf5ea8SMatthew G. Knepley   for (i = 0; i < patch->npatch; ++i) {
15644bbf5ea8SMatthew G. Knepley     ierr = KSPSetUp(patch->ksp[i]);CHKERRQ(ierr);
15654bbf5ea8SMatthew G. Knepley     ierr = KSPGetConvergedReason(patch->ksp[i], &reason);CHKERRQ(ierr);
15664bbf5ea8SMatthew G. Knepley     if (reason == KSP_DIVERGED_PCSETUP_FAILED) pc->failedreason = PC_SUBPC_ERROR;
15674bbf5ea8SMatthew G. Knepley   }
15684bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
15694bbf5ea8SMatthew G. Knepley }
15704bbf5ea8SMatthew G. Knepley 
15714bbf5ea8SMatthew G. Knepley static PetscErrorCode PCView_PATCH(PC pc, PetscViewer viewer)
15724bbf5ea8SMatthew G. Knepley {
15734bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch = (PC_PATCH *) pc->data;
15744bbf5ea8SMatthew G. Knepley   PetscViewer    sviewer;
15754bbf5ea8SMatthew G. Knepley   PetscBool      isascii;
15764bbf5ea8SMatthew G. Knepley   PetscMPIInt    rank;
15774bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
15784bbf5ea8SMatthew G. Knepley 
15794bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
15804bbf5ea8SMatthew G. Knepley   /* TODO Redo tabbing with set tbas in new style */
15814bbf5ea8SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
15824bbf5ea8SMatthew G. Knepley   if (!isascii) PetscFunctionReturn(0);
15834bbf5ea8SMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) pc), &rank);CHKERRQ(ierr);
15844bbf5ea8SMatthew G. Knepley   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15854bbf5ea8SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Subspace Correction preconditioner with %d patches\n", patch->npatch);CHKERRQ(ierr);
15864bbf5ea8SMatthew G. Knepley   if (patch->multiplicative) {ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: multiplicative\n");CHKERRQ(ierr);}
15874bbf5ea8SMatthew G. Knepley   else                       {ierr = PetscViewerASCIIPrintf(viewer, "Schwarz type: additive\n");CHKERRQ(ierr);}
15884bbf5ea8SMatthew G. Knepley   if (patch->partition_of_unity) {ierr = PetscViewerASCIIPrintf(viewer, "Weighting by partition of unity\n");CHKERRQ(ierr);}
15894bbf5ea8SMatthew G. Knepley   else                           {ierr = PetscViewerASCIIPrintf(viewer, "Not weighting by partition of unity\n");CHKERRQ(ierr);}
15904bbf5ea8SMatthew G. Knepley   if (patch->symmetrise_sweep) {ierr = PetscViewerASCIIPrintf(viewer, "Symmetrising sweep (start->end, then end->start)\n");CHKERRQ(ierr);}
15914bbf5ea8SMatthew G. Knepley   else                         {ierr = PetscViewerASCIIPrintf(viewer, "Not symmetrising sweep\n");CHKERRQ(ierr);}
15924bbf5ea8SMatthew G. Knepley   if (!patch->save_operators) {ierr = PetscViewerASCIIPrintf(viewer, "Not saving patch operators (rebuilt every PCApply)\n");CHKERRQ(ierr);}
15934bbf5ea8SMatthew G. Knepley   else                        {ierr = PetscViewerASCIIPrintf(viewer, "Saving patch operators (rebuilt every PCSetUp)\n");CHKERRQ(ierr);}
15944bbf5ea8SMatthew G. Knepley   if (patch->patchconstructop == PCPatchConstruct_Star)       {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: star\n");CHKERRQ(ierr);}
15954bbf5ea8SMatthew G. Knepley   else if (patch->patchconstructop == PCPatchConstruct_Vanka) {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: Vanka\n");CHKERRQ(ierr);}
15964bbf5ea8SMatthew G. Knepley   else if (patch->patchconstructop == PCPatchConstruct_User)  {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: user-specified\n");CHKERRQ(ierr);}
15974bbf5ea8SMatthew G. Knepley   else                                                        {ierr = PetscViewerASCIIPrintf(viewer, "Patch construction operator: unknown\n");CHKERRQ(ierr);}
15984bbf5ea8SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "KSP on patches (all same):\n");CHKERRQ(ierr);
15994bbf5ea8SMatthew G. Knepley   if (patch->ksp) {
16004bbf5ea8SMatthew G. Knepley     ierr = PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
16014bbf5ea8SMatthew G. Knepley     if (!rank) {
16024bbf5ea8SMatthew G. Knepley       ierr = PetscViewerASCIIPushTab(sviewer);CHKERRQ(ierr);
16034bbf5ea8SMatthew G. Knepley       ierr = KSPView(patch->ksp[0], sviewer);CHKERRQ(ierr);
16044bbf5ea8SMatthew G. Knepley       ierr = PetscViewerASCIIPopTab(sviewer);CHKERRQ(ierr);
16054bbf5ea8SMatthew G. Knepley     }
16064bbf5ea8SMatthew G. Knepley     ierr = PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer);CHKERRQ(ierr);
16074bbf5ea8SMatthew G. Knepley   } else {
16084bbf5ea8SMatthew G. Knepley     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
16094bbf5ea8SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "KSP not yet set.\n");CHKERRQ(ierr);
16104bbf5ea8SMatthew G. Knepley     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
16114bbf5ea8SMatthew G. Knepley   }
16124bbf5ea8SMatthew G. Knepley   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
16134bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
16144bbf5ea8SMatthew G. Knepley }
16154bbf5ea8SMatthew G. Knepley 
1616*642283e9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCCreate_Patch(PC pc)
16174bbf5ea8SMatthew G. Knepley {
16184bbf5ea8SMatthew G. Knepley   PC_PATCH      *patch;
16194bbf5ea8SMatthew G. Knepley   PetscErrorCode ierr;
16204bbf5ea8SMatthew G. Knepley 
16214bbf5ea8SMatthew G. Knepley   PetscFunctionBegin;
16224bbf5ea8SMatthew G. Knepley   ierr = PetscNewLog(pc, &patch);CHKERRQ(ierr);
16234bbf5ea8SMatthew G. Knepley 
16244bbf5ea8SMatthew G. Knepley   /* Set some defaults */
16254bbf5ea8SMatthew G. Knepley   patch->sub_mat_type       = NULL;
16264bbf5ea8SMatthew G. Knepley   patch->save_operators     = PETSC_TRUE;
16274bbf5ea8SMatthew G. Knepley   patch->partition_of_unity = PETSC_FALSE;
16284bbf5ea8SMatthew G. Knepley   patch->multiplicative     = PETSC_FALSE;
16294bbf5ea8SMatthew G. Knepley   patch->codim              = -1;
16304bbf5ea8SMatthew G. Knepley   patch->dim                = -1;
16314bbf5ea8SMatthew G. Knepley   patch->exclude_subspace   = -1;
16324bbf5ea8SMatthew G. Knepley   patch->vankadim           = -1;
16334bbf5ea8SMatthew G. Knepley   patch->patchconstructop   = PCPatchConstruct_Star;
16344bbf5ea8SMatthew G. Knepley   patch->print_patches      = PETSC_FALSE;
16354bbf5ea8SMatthew G. Knepley   patch->symmetrise_sweep   = PETSC_FALSE;
16364bbf5ea8SMatthew G. Knepley   patch->nuserIS            = 0;
16374bbf5ea8SMatthew G. Knepley   patch->userIS             = NULL;
16384bbf5ea8SMatthew G. Knepley   patch->iterationSet       = NULL;
16394bbf5ea8SMatthew G. Knepley   patch->user_patches       = PETSC_FALSE;
16404bbf5ea8SMatthew G. Knepley 
16414bbf5ea8SMatthew G. Knepley   pc->data                 = (void *) patch;
16424bbf5ea8SMatthew G. Knepley   pc->ops->apply           = PCApply_PATCH;
16434bbf5ea8SMatthew G. Knepley   pc->ops->applytranspose  = 0; /* PCApplyTranspose_PATCH; */
16444bbf5ea8SMatthew G. Knepley   pc->ops->setup           = PCSetUp_PATCH;
16454bbf5ea8SMatthew G. Knepley   pc->ops->reset           = PCReset_PATCH;
16464bbf5ea8SMatthew G. Knepley   pc->ops->destroy         = PCDestroy_PATCH;
16474bbf5ea8SMatthew G. Knepley   pc->ops->setfromoptions  = PCSetFromOptions_PATCH;
16484bbf5ea8SMatthew G. Knepley   pc->ops->setuponblocks   = PCSetUpOnBlocks_PATCH;
16494bbf5ea8SMatthew G. Knepley   pc->ops->view            = PCView_PATCH;
16504bbf5ea8SMatthew G. Knepley   pc->ops->applyrichardson = 0;
16514bbf5ea8SMatthew G. Knepley 
16524bbf5ea8SMatthew G. Knepley   PetscFunctionReturn(0);
16534bbf5ea8SMatthew G. Knepley }
1654