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