1*6858538eSMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 2*6858538eSMatthew G. Knepley 3*6858538eSMatthew G. Knepley #include <petscdmplex.h> 4*6858538eSMatthew G. Knepley 5*6858538eSMatthew G. Knepley /*@C 6*6858538eSMatthew G. Knepley DMGetPeriodicity - Get the description of mesh periodicity 7*6858538eSMatthew G. Knepley 8*6858538eSMatthew G. Knepley Input Parameter: 9*6858538eSMatthew G. Knepley . dm - The DM object 10*6858538eSMatthew G. Knepley 11*6858538eSMatthew G. Knepley Output Parameters: 12*6858538eSMatthew G. Knepley + maxCell - Over distances greater than this, we can assume a point has crossed over to another sheet, when trying to localize cell coordinates 13*6858538eSMatthew G. Knepley - L - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0 14*6858538eSMatthew G. Knepley 15*6858538eSMatthew G. Knepley Level: developer 16*6858538eSMatthew G. Knepley 17*6858538eSMatthew G. Knepley .seealso: `DMGetPeriodicity()` 18*6858538eSMatthew G. Knepley @*/ 19*6858538eSMatthew G. Knepley PetscErrorCode DMGetPeriodicity(DM dm, const PetscReal **maxCell, const PetscReal **L) 20*6858538eSMatthew G. Knepley { 21*6858538eSMatthew G. Knepley PetscFunctionBegin; 22*6858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 23*6858538eSMatthew G. Knepley if (L) *L = dm->L; 24*6858538eSMatthew G. Knepley if (maxCell) *maxCell = dm->maxCell; 25*6858538eSMatthew G. Knepley PetscFunctionReturn(0); 26*6858538eSMatthew G. Knepley } 27*6858538eSMatthew G. Knepley 28*6858538eSMatthew G. Knepley /*@C 29*6858538eSMatthew G. Knepley DMSetPeriodicity - Set the description of mesh periodicity 30*6858538eSMatthew G. Knepley 31*6858538eSMatthew G. Knepley Input Parameters: 32*6858538eSMatthew G. Knepley + dm - The DM object 33*6858538eSMatthew G. Knepley . maxCell - Over distances greater than this, we can assume a point has crossed over to another sheet, when trying to localize cell coordinates. Pass NULL to remove such information. 34*6858538eSMatthew G. Knepley - L - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0 35*6858538eSMatthew G. Knepley 36*6858538eSMatthew G. Knepley Level: developer 37*6858538eSMatthew G. Knepley 38*6858538eSMatthew G. Knepley .seealso: `DMGetPeriodicity()` 39*6858538eSMatthew G. Knepley @*/ 40*6858538eSMatthew G. Knepley PetscErrorCode DMSetPeriodicity(DM dm, const PetscReal maxCell[], const PetscReal L[]) 41*6858538eSMatthew G. Knepley { 42*6858538eSMatthew G. Knepley PetscInt dim, d; 43*6858538eSMatthew G. Knepley 44*6858538eSMatthew G. Knepley PetscFunctionBegin; 45*6858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm,DM_CLASSID,1); 46*6858538eSMatthew G. Knepley if (maxCell) {PetscValidRealPointer(maxCell,2);} 47*6858538eSMatthew G. Knepley if (L) {PetscValidRealPointer(L,3);} 48*6858538eSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 49*6858538eSMatthew G. Knepley if (maxCell) { 50*6858538eSMatthew G. Knepley if (!dm->maxCell) PetscCall(PetscMalloc1(dim, &dm->maxCell)); 51*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) dm->maxCell[d] = maxCell[d]; 52*6858538eSMatthew G. Knepley } else { /* remove maxCell information to disable automatic computation of localized vertices */ 53*6858538eSMatthew G. Knepley PetscCall(PetscFree(dm->maxCell)); 54*6858538eSMatthew G. Knepley dm->maxCell = NULL; 55*6858538eSMatthew G. Knepley } 56*6858538eSMatthew G. Knepley if (L) { 57*6858538eSMatthew G. Knepley if (!dm->L) PetscCall(PetscMalloc1(dim, &dm->L)); 58*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) dm->L[d] = L[d]; 59*6858538eSMatthew G. Knepley } else { /* remove L information to disable automatic computation of localized vertices */ 60*6858538eSMatthew G. Knepley PetscCall(PetscFree(dm->L)); 61*6858538eSMatthew G. Knepley dm->L = NULL; 62*6858538eSMatthew G. Knepley } 63*6858538eSMatthew G. Knepley PetscCheck((dm->maxCell && dm->L) || (!dm->maxCell && !dm->L), PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Cannot set only one of maxCell/L"); 64*6858538eSMatthew G. Knepley PetscFunctionReturn(0); 65*6858538eSMatthew G. Knepley } 66*6858538eSMatthew G. Knepley 67*6858538eSMatthew G. Knepley /*@ 68*6858538eSMatthew G. Knepley DMLocalizeCoordinate - If a mesh is periodic (a torus with lengths L_i, some of which can be infinite), project the coordinate onto [0, L_i) in each dimension. 69*6858538eSMatthew G. Knepley 70*6858538eSMatthew G. Knepley Input Parameters: 71*6858538eSMatthew G. Knepley + dm - The DM 72*6858538eSMatthew G. Knepley . in - The input coordinate point (dim numbers) 73*6858538eSMatthew G. Knepley - endpoint - Include the endpoint L_i 74*6858538eSMatthew G. Knepley 75*6858538eSMatthew G. Knepley Output Parameter: 76*6858538eSMatthew G. Knepley . out - The localized coordinate point 77*6858538eSMatthew G. Knepley 78*6858538eSMatthew G. Knepley Level: developer 79*6858538eSMatthew G. Knepley 80*6858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()` 81*6858538eSMatthew G. Knepley @*/ 82*6858538eSMatthew G. Knepley PetscErrorCode DMLocalizeCoordinate(DM dm, const PetscScalar in[], PetscBool endpoint, PetscScalar out[]) 83*6858538eSMatthew G. Knepley { 84*6858538eSMatthew G. Knepley PetscInt dim, d; 85*6858538eSMatthew G. Knepley 86*6858538eSMatthew G. Knepley PetscFunctionBegin; 87*6858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &dim)); 88*6858538eSMatthew G. Knepley if (!dm->maxCell) { 89*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) out[d] = in[d]; 90*6858538eSMatthew G. Knepley } else { 91*6858538eSMatthew G. Knepley if (endpoint) { 92*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 93*6858538eSMatthew G. Knepley if ((PetscAbsReal(PetscRealPart(in[d])/dm->L[d] - PetscFloorReal(PetscRealPart(in[d])/dm->L[d])) < PETSC_SMALL) && (PetscRealPart(in[d])/dm->L[d] > PETSC_SMALL)) { 94*6858538eSMatthew G. Knepley out[d] = in[d] - dm->L[d]*(PetscFloorReal(PetscRealPart(in[d])/dm->L[d]) - 1); 95*6858538eSMatthew G. Knepley } else { 96*6858538eSMatthew G. Knepley out[d] = in[d] - dm->L[d]*PetscFloorReal(PetscRealPart(in[d])/dm->L[d]); 97*6858538eSMatthew G. Knepley } 98*6858538eSMatthew G. Knepley } 99*6858538eSMatthew G. Knepley } else { 100*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 101*6858538eSMatthew G. Knepley out[d] = in[d] - dm->L[d]*PetscFloorReal(PetscRealPart(in[d])/dm->L[d]); 102*6858538eSMatthew G. Knepley } 103*6858538eSMatthew G. Knepley } 104*6858538eSMatthew G. Knepley } 105*6858538eSMatthew G. Knepley PetscFunctionReturn(0); 106*6858538eSMatthew G. Knepley } 107*6858538eSMatthew G. Knepley 108*6858538eSMatthew G. Knepley /* 109*6858538eSMatthew G. Knepley DMLocalizeCoordinate_Internal - If a mesh is periodic, and the input point is far from the anchor, pick the coordinate sheet of the torus which moves it closer. 110*6858538eSMatthew G. Knepley 111*6858538eSMatthew G. Knepley Input Parameters: 112*6858538eSMatthew G. Knepley + dm - The DM 113*6858538eSMatthew G. Knepley . dim - The spatial dimension 114*6858538eSMatthew G. Knepley . anchor - The anchor point, the input point can be no more than maxCell away from it 115*6858538eSMatthew G. Knepley - in - The input coordinate point (dim numbers) 116*6858538eSMatthew G. Knepley 117*6858538eSMatthew G. Knepley Output Parameter: 118*6858538eSMatthew G. Knepley . out - The localized coordinate point 119*6858538eSMatthew G. Knepley 120*6858538eSMatthew G. Knepley Level: developer 121*6858538eSMatthew G. Knepley 122*6858538eSMatthew G. Knepley Note: This is meant to get a set of coordinates close to each other, as in a cell. The anchor is usually the one of the vertices on a containing cell 123*6858538eSMatthew G. Knepley 124*6858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()` 125*6858538eSMatthew G. Knepley */ 126*6858538eSMatthew G. Knepley PetscErrorCode DMLocalizeCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[]) 127*6858538eSMatthew G. Knepley { 128*6858538eSMatthew G. Knepley PetscInt d; 129*6858538eSMatthew G. Knepley 130*6858538eSMatthew G. Knepley PetscFunctionBegin; 131*6858538eSMatthew G. Knepley if (!dm->maxCell) { 132*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) out[d] = in[d]; 133*6858538eSMatthew G. Knepley } else { 134*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 135*6858538eSMatthew G. Knepley if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > dm->maxCell[d])) { 136*6858538eSMatthew G. Knepley out[d] = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d]; 137*6858538eSMatthew G. Knepley } else { 138*6858538eSMatthew G. Knepley out[d] = in[d]; 139*6858538eSMatthew G. Knepley } 140*6858538eSMatthew G. Knepley } 141*6858538eSMatthew G. Knepley } 142*6858538eSMatthew G. Knepley PetscFunctionReturn(0); 143*6858538eSMatthew G. Knepley } 144*6858538eSMatthew G. Knepley 145*6858538eSMatthew G. Knepley PetscErrorCode DMLocalizeCoordinateReal_Internal(DM dm, PetscInt dim, const PetscReal anchor[], const PetscReal in[], PetscReal out[]) 146*6858538eSMatthew G. Knepley { 147*6858538eSMatthew G. Knepley PetscInt d; 148*6858538eSMatthew G. Knepley 149*6858538eSMatthew G. Knepley PetscFunctionBegin; 150*6858538eSMatthew G. Knepley if (!dm->maxCell) { 151*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) out[d] = in[d]; 152*6858538eSMatthew G. Knepley } else { 153*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 154*6858538eSMatthew G. Knepley if ((dm->L[d] > 0.0) && (PetscAbsReal(anchor[d] - in[d]) > dm->maxCell[d])) { 155*6858538eSMatthew G. Knepley out[d] = anchor[d] > in[d] ? dm->L[d] + in[d] : in[d] - dm->L[d]; 156*6858538eSMatthew G. Knepley } else { 157*6858538eSMatthew G. Knepley out[d] = in[d]; 158*6858538eSMatthew G. Knepley } 159*6858538eSMatthew G. Knepley } 160*6858538eSMatthew G. Knepley } 161*6858538eSMatthew G. Knepley PetscFunctionReturn(0); 162*6858538eSMatthew G. Knepley } 163*6858538eSMatthew G. Knepley 164*6858538eSMatthew G. Knepley /* 165*6858538eSMatthew G. Knepley DMLocalizeAddCoordinate_Internal - If a mesh is periodic, and the input point is far from the anchor, pick the coordinate sheet of the torus which moves it closer. 166*6858538eSMatthew G. Knepley 167*6858538eSMatthew G. Knepley Input Parameters: 168*6858538eSMatthew G. Knepley + dm - The DM 169*6858538eSMatthew G. Knepley . dim - The spatial dimension 170*6858538eSMatthew G. Knepley . anchor - The anchor point, the input point can be no more than maxCell away from it 171*6858538eSMatthew G. Knepley . in - The input coordinate delta (dim numbers) 172*6858538eSMatthew G. Knepley - out - The input coordinate point (dim numbers) 173*6858538eSMatthew G. Knepley 174*6858538eSMatthew G. Knepley Output Parameter: 175*6858538eSMatthew G. Knepley . out - The localized coordinate in + out 176*6858538eSMatthew G. Knepley 177*6858538eSMatthew G. Knepley Level: developer 178*6858538eSMatthew G. Knepley 179*6858538eSMatthew G. Knepley Note: This is meant to get a set of coordinates close to each other, as in a cell. The anchor is usually the one of the vertices on a containing cell 180*6858538eSMatthew G. Knepley 181*6858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMLocalizeCoordinate()` 182*6858538eSMatthew G. Knepley */ 183*6858538eSMatthew G. Knepley PetscErrorCode DMLocalizeAddCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[]) 184*6858538eSMatthew G. Knepley { 185*6858538eSMatthew G. Knepley PetscInt d; 186*6858538eSMatthew G. Knepley 187*6858538eSMatthew G. Knepley PetscFunctionBegin; 188*6858538eSMatthew G. Knepley if (!dm->maxCell) { 189*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) out[d] += in[d]; 190*6858538eSMatthew G. Knepley } else { 191*6858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 192*6858538eSMatthew G. Knepley const PetscReal maxC = dm->maxCell[d]; 193*6858538eSMatthew G. Knepley 194*6858538eSMatthew G. Knepley if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > maxC)) { 195*6858538eSMatthew G. Knepley const PetscScalar newCoord = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d]; 196*6858538eSMatthew G. Knepley 197*6858538eSMatthew G. Knepley if (PetscAbsScalar(newCoord - anchor[d]) > maxC) 198*6858538eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT "-Coordinate %g more than %g away from anchor %g", d, (double) PetscRealPart(in[d]), (double) maxC, (double) PetscRealPart(anchor[d])); 199*6858538eSMatthew G. Knepley out[d] += newCoord; 200*6858538eSMatthew G. Knepley } else { 201*6858538eSMatthew G. Knepley out[d] += in[d]; 202*6858538eSMatthew G. Knepley } 203*6858538eSMatthew G. Knepley } 204*6858538eSMatthew G. Knepley } 205*6858538eSMatthew G. Knepley PetscFunctionReturn(0); 206*6858538eSMatthew G. Knepley } 207*6858538eSMatthew G. Knepley 208*6858538eSMatthew G. Knepley /*@ 209*6858538eSMatthew G. Knepley DMGetCoordinatesLocalizedLocal - Check if the DM coordinates have been localized for cells on this process 210*6858538eSMatthew G. Knepley 211*6858538eSMatthew G. Knepley Not collective 212*6858538eSMatthew G. Knepley 213*6858538eSMatthew G. Knepley Input Parameter: 214*6858538eSMatthew G. Knepley . dm - The DM 215*6858538eSMatthew G. Knepley 216*6858538eSMatthew G. Knepley Output Parameter: 217*6858538eSMatthew G. Knepley areLocalized - True if localized 218*6858538eSMatthew G. Knepley 219*6858538eSMatthew G. Knepley Level: developer 220*6858538eSMatthew G. Knepley 221*6858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMGetCoordinatesLocalized()`, `DMSetPeriodicity()` 222*6858538eSMatthew G. Knepley @*/ 223*6858538eSMatthew G. Knepley PetscErrorCode DMGetCoordinatesLocalizedLocal(DM dm, PetscBool *areLocalized) 224*6858538eSMatthew G. Knepley { 225*6858538eSMatthew G. Knepley PetscFunctionBegin; 226*6858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 227*6858538eSMatthew G. Knepley PetscValidBoolPointer(areLocalized, 2); 228*6858538eSMatthew G. Knepley *areLocalized = dm->coordinates[1].dim < 0 ? PETSC_FALSE : PETSC_TRUE; 229*6858538eSMatthew G. Knepley PetscFunctionReturn(0); 230*6858538eSMatthew G. Knepley } 231*6858538eSMatthew G. Knepley 232*6858538eSMatthew G. Knepley /*@ 233*6858538eSMatthew G. Knepley DMGetCoordinatesLocalized - Check if the DM coordinates have been localized for cells 234*6858538eSMatthew G. Knepley 235*6858538eSMatthew G. Knepley Collective on dm 236*6858538eSMatthew G. Knepley 237*6858538eSMatthew G. Knepley Input Parameter: 238*6858538eSMatthew G. Knepley . dm - The DM 239*6858538eSMatthew G. Knepley 240*6858538eSMatthew G. Knepley Output Parameter: 241*6858538eSMatthew G. Knepley areLocalized - True if localized 242*6858538eSMatthew G. Knepley 243*6858538eSMatthew G. Knepley Level: developer 244*6858538eSMatthew G. Knepley 245*6858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMSetPeriodicity()`, `DMGetCoordinatesLocalizedLocal()` 246*6858538eSMatthew G. Knepley @*/ 247*6858538eSMatthew G. Knepley PetscErrorCode DMGetCoordinatesLocalized(DM dm,PetscBool *areLocalized) 248*6858538eSMatthew G. Knepley { 249*6858538eSMatthew G. Knepley PetscBool localized; 250*6858538eSMatthew G. Knepley 251*6858538eSMatthew G. Knepley PetscFunctionBegin; 252*6858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 253*6858538eSMatthew G. Knepley PetscValidBoolPointer(areLocalized, 2); 254*6858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalizedLocal(dm, &localized)); 255*6858538eSMatthew G. Knepley PetscCall(MPIU_Allreduce(&localized, areLocalized, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject) dm))); 256*6858538eSMatthew G. Knepley PetscFunctionReturn(0); 257*6858538eSMatthew G. Knepley } 258*6858538eSMatthew G. Knepley 259*6858538eSMatthew G. Knepley /*@ 260*6858538eSMatthew G. Knepley DMLocalizeCoordinates - If a mesh is periodic, create local coordinates for cells having periodic faces 261*6858538eSMatthew G. Knepley 262*6858538eSMatthew G. Knepley Collective on dm 263*6858538eSMatthew G. Knepley 264*6858538eSMatthew G. Knepley Input Parameter: 265*6858538eSMatthew G. Knepley . dm - The DM 266*6858538eSMatthew G. Knepley 267*6858538eSMatthew G. Knepley Level: developer 268*6858538eSMatthew G. Knepley 269*6858538eSMatthew G. Knepley .seealso: `DMSetPeriodicity()`, `DMLocalizeCoordinate()`, `DMLocalizeAddCoordinate()` 270*6858538eSMatthew G. Knepley @*/ 271*6858538eSMatthew G. Knepley PetscErrorCode DMLocalizeCoordinates(DM dm) 272*6858538eSMatthew G. Knepley { 273*6858538eSMatthew G. Knepley DM cdm, cdgdm, cplex, plex; 274*6858538eSMatthew G. Knepley PetscSection cs, csDG; 275*6858538eSMatthew G. Knepley Vec coordinates, cVec; 276*6858538eSMatthew G. Knepley PetscScalar *coordsDG, *anchor, *localized; 277*6858538eSMatthew G. Knepley const PetscReal *L; 278*6858538eSMatthew G. Knepley PetscInt Nc, vStart, vEnd, sStart, sEnd, newStart = PETSC_MAX_INT, newEnd = PETSC_MIN_INT, bs, coordSize; 279*6858538eSMatthew G. Knepley PetscBool isLocalized, sparseLocalize = dm->sparseLocalize, useDG = PETSC_FALSE, useDGGlobal; 280*6858538eSMatthew G. Knepley PetscInt maxHeight = 0, h; 281*6858538eSMatthew G. Knepley PetscInt *pStart = NULL, *pEnd = NULL; 282*6858538eSMatthew G. Knepley MPI_Comm comm; 283*6858538eSMatthew G. Knepley 284*6858538eSMatthew G. Knepley PetscFunctionBegin; 285*6858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 286*6858538eSMatthew G. Knepley PetscCall(DMGetPeriodicity(dm, NULL, &L)); 287*6858538eSMatthew G. Knepley /* Cannot automatically localize without L and maxCell right now */ 288*6858538eSMatthew G. Knepley if (!L) PetscFunctionReturn(0); 289*6858538eSMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 290*6858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalized(dm, &isLocalized)); 291*6858538eSMatthew G. Knepley if (isLocalized) PetscFunctionReturn(0); 292*6858538eSMatthew G. Knepley 293*6858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 294*6858538eSMatthew G. Knepley PetscCall(DMConvert(dm, DMPLEX, &plex)); 295*6858538eSMatthew G. Knepley PetscCall(DMConvert(cdm, DMPLEX, &cplex)); 296*6858538eSMatthew G. Knepley if (cplex) { 297*6858538eSMatthew G. Knepley PetscCall(DMPlexGetDepthStratum(cplex, 0, &vStart, &vEnd)); 298*6858538eSMatthew G. Knepley PetscCall(DMPlexGetMaxProjectionHeight(cplex, &maxHeight)); 299*6858538eSMatthew G. Knepley PetscCall(DMGetWorkArray(dm, 2*(maxHeight + 1), MPIU_INT, &pStart)); 300*6858538eSMatthew G. Knepley pEnd = &pStart[maxHeight + 1]; 301*6858538eSMatthew G. Knepley newStart = vStart; 302*6858538eSMatthew G. Knepley newEnd = vEnd; 303*6858538eSMatthew G. Knepley for (h = 0; h <= maxHeight; h++) { 304*6858538eSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(cplex, h, &pStart[h], &pEnd[h])); 305*6858538eSMatthew G. Knepley newStart = PetscMin(newStart, pStart[h]); 306*6858538eSMatthew G. Knepley newEnd = PetscMax(newEnd, pEnd[h]); 307*6858538eSMatthew G. Knepley } 308*6858538eSMatthew G. Knepley } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Coordinate localization requires a DMPLEX coordinate DM"); 309*6858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 310*6858538eSMatthew G. Knepley PetscCheck(coordinates, comm, PETSC_ERR_SUP, "Missing local coordinate vector"); 311*6858538eSMatthew G. Knepley PetscCall(DMGetCoordinateSection(dm, &cs)); 312*6858538eSMatthew G. Knepley PetscCall(VecGetBlockSize(coordinates, &bs)); 313*6858538eSMatthew G. Knepley PetscCall(PetscSectionGetChart(cs, &sStart, &sEnd)); 314*6858538eSMatthew G. Knepley 315*6858538eSMatthew G. Knepley PetscCall(PetscSectionCreate(comm, &csDG)); 316*6858538eSMatthew G. Knepley PetscCall(PetscSectionSetNumFields(csDG, 1)); 317*6858538eSMatthew G. Knepley PetscCall(PetscSectionGetFieldComponents(cs, 0, &Nc)); 318*6858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldComponents(csDG, 0, Nc)); 319*6858538eSMatthew G. Knepley PetscCall(PetscSectionSetChart(csDG, newStart, newEnd)); 320*6858538eSMatthew G. Knepley PetscCheck(bs == Nc, comm, PETSC_ERR_ARG_INCOMP, "Coordinate block size %" PetscInt_FMT " != %" PetscInt_FMT " number of components", bs, Nc); 321*6858538eSMatthew G. Knepley 322*6858538eSMatthew G. Knepley PetscCall(DMGetWorkArray(dm, 2 * Nc, MPIU_SCALAR, &anchor)); 323*6858538eSMatthew G. Knepley localized = &anchor[Nc]; 324*6858538eSMatthew G. Knepley for (h = 0; h <= maxHeight; h++) { 325*6858538eSMatthew G. Knepley PetscInt cStart = pStart[h], cEnd = pEnd[h], c; 326*6858538eSMatthew G. Knepley 327*6858538eSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 328*6858538eSMatthew G. Knepley PetscScalar *cellCoords = NULL; 329*6858538eSMatthew G. Knepley DMPolytopeType ct; 330*6858538eSMatthew G. Knepley PetscInt dof, d, p; 331*6858538eSMatthew G. Knepley 332*6858538eSMatthew G. Knepley PetscCall(DMPlexGetCellType(plex, c, &ct)); 333*6858538eSMatthew G. Knepley if (ct == DM_POLYTOPE_FV_GHOST) continue; 334*6858538eSMatthew G. Knepley PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords)); 335*6858538eSMatthew G. Knepley PetscCheck(!(dof % Nc), comm, PETSC_ERR_ARG_INCOMP, "Coordinate size on cell %" PetscInt_FMT " closure %" PetscInt_FMT " not divisible by %" PetscInt_FMT " number of components", c, dof, Nc); 336*6858538eSMatthew G. Knepley for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[d]; 337*6858538eSMatthew G. Knepley for (p = 0; p < dof/Nc; ++p) { 338*6858538eSMatthew G. Knepley PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p*Nc], localized)); 339*6858538eSMatthew G. Knepley for (d = 0; d < Nc; ++d) if (cellCoords[p*Nc + d] != localized[d]) break; 340*6858538eSMatthew G. Knepley if (d < Nc) break; 341*6858538eSMatthew G. Knepley } 342*6858538eSMatthew G. Knepley if (p < dof/Nc) useDG = PETSC_TRUE; 343*6858538eSMatthew G. Knepley if (p < dof/Nc || !sparseLocalize) { 344*6858538eSMatthew G. Knepley PetscCall(PetscSectionSetDof(csDG, c, dof)); 345*6858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(csDG, c, 0, dof)); 346*6858538eSMatthew G. Knepley } 347*6858538eSMatthew G. Knepley PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords)); 348*6858538eSMatthew G. Knepley } 349*6858538eSMatthew G. Knepley } 350*6858538eSMatthew G. Knepley PetscCallMPI(MPI_Allreduce(&useDG, &useDGGlobal, 1, MPIU_BOOL, MPI_LOR, comm)); 351*6858538eSMatthew G. Knepley if (!useDGGlobal) goto end; 352*6858538eSMatthew G. Knepley 353*6858538eSMatthew G. Knepley PetscCall(PetscSectionSetUp(csDG)); 354*6858538eSMatthew G. Knepley PetscCall(PetscSectionGetStorageSize(csDG, &coordSize)); 355*6858538eSMatthew G. Knepley PetscCall(VecCreate(PETSC_COMM_SELF, &cVec)); 356*6858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject) cVec, "coordinates")); 357*6858538eSMatthew G. Knepley PetscCall(VecSetBlockSize(cVec, bs)); 358*6858538eSMatthew G. Knepley PetscCall(VecSetSizes(cVec, coordSize, PETSC_DETERMINE)); 359*6858538eSMatthew G. Knepley PetscCall(VecSetType(cVec, VECSTANDARD)); 360*6858538eSMatthew G. Knepley PetscCall(VecGetArray(cVec, &coordsDG)); 361*6858538eSMatthew G. Knepley for (h = 0; h <= maxHeight; h++) { 362*6858538eSMatthew G. Knepley PetscInt cStart = pStart[h], cEnd = pEnd[h], c; 363*6858538eSMatthew G. Knepley 364*6858538eSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 365*6858538eSMatthew G. Knepley PetscScalar *cellCoords = NULL; 366*6858538eSMatthew G. Knepley PetscInt p = 0, q, dof, cdof, d, offDG; 367*6858538eSMatthew G. Knepley 368*6858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(csDG, c, &cdof)); 369*6858538eSMatthew G. Knepley if (!cdof) continue; 370*6858538eSMatthew G. Knepley PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords)); 371*6858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(csDG, c, &offDG)); 372*6858538eSMatthew G. Knepley // We need the cell to fit into [0, L] 373*6858538eSMatthew G. Knepley for (q = 0; q < dof/Nc; ++q) { 374*6858538eSMatthew G. Knepley // Select a trial anchor 375*6858538eSMatthew G. Knepley for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[q*Nc+d]; 376*6858538eSMatthew G. Knepley for (p = 0; p < dof/Nc; ++p) { 377*6858538eSMatthew G. Knepley PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p*Nc], &coordsDG[offDG + p*Nc])); 378*6858538eSMatthew G. Knepley for (d = 0; d < Nc; ++d) 379*6858538eSMatthew G. Knepley if (L[d] > 0. && ((PetscRealPart(coordsDG[offDG + p*Nc + d]) < 0.) || (PetscRealPart(coordsDG[offDG + p*Nc + d]) > L[d]))) break; 380*6858538eSMatthew G. Knepley if (d < Nc) break; 381*6858538eSMatthew G. Knepley } 382*6858538eSMatthew G. Knepley if (p == dof/Nc) break; 383*6858538eSMatthew G. Knepley } 384*6858538eSMatthew G. Knepley PetscCheck(p == dof/Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " does not fit into the torus [0, L]", c); 385*6858538eSMatthew G. Knepley PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords)); 386*6858538eSMatthew G. Knepley } 387*6858538eSMatthew G. Knepley } 388*6858538eSMatthew G. Knepley PetscCall(VecRestoreArray(cVec, &coordsDG)); 389*6858538eSMatthew G. Knepley PetscCall(DMClone(cdm, &cdgdm)); 390*6858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinateDM(dm, cdgdm)); 391*6858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinateSection(dm, PETSC_DETERMINE, csDG)); 392*6858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinatesLocal(dm, cVec)); 393*6858538eSMatthew G. Knepley PetscCall(VecDestroy(&cVec)); 394*6858538eSMatthew G. Knepley PetscCall(DMDestroy(&cdgdm)); 395*6858538eSMatthew G. Knepley 396*6858538eSMatthew G. Knepley end: 397*6858538eSMatthew G. Knepley PetscCall(DMRestoreWorkArray(dm, 2 * bs, MPIU_SCALAR, &anchor)); 398*6858538eSMatthew G. Knepley PetscCall(DMRestoreWorkArray(dm, 2*(maxHeight + 1), MPIU_INT, &pStart)); 399*6858538eSMatthew G. Knepley PetscCall(PetscSectionDestroy(&csDG)); 400*6858538eSMatthew G. Knepley PetscCall(DMDestroy(&plex)); 401*6858538eSMatthew G. Knepley PetscCall(DMDestroy(&cplex)); 402*6858538eSMatthew G. Knepley PetscFunctionReturn(0); 403*6858538eSMatthew G. Knepley } 404