120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac #include <petscdmplex.h> 320cf1dd8SToby Isaac 4*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp) 520cf1dd8SToby Isaac { 620cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 7*6f905325SMatthew G. Knepley PetscInt i; 8*6f905325SMatthew G. Knepley PetscErrorCode ierr; 920cf1dd8SToby Isaac 1020cf1dd8SToby Isaac PetscFunctionBegin; 11*6f905325SMatthew G. Knepley if (lag->symmetries) { 12*6f905325SMatthew G. Knepley PetscInt **selfSyms = lag->symmetries[0]; 13*6f905325SMatthew G. Knepley 14*6f905325SMatthew G. Knepley if (selfSyms) { 15*6f905325SMatthew G. Knepley PetscInt i, **allocated = &selfSyms[-lag->selfSymOff]; 16*6f905325SMatthew G. Knepley 17*6f905325SMatthew G. Knepley for (i = 0; i < lag->numSelfSym; i++) { 18*6f905325SMatthew G. Knepley ierr = PetscFree(allocated[i]);CHKERRQ(ierr); 19*6f905325SMatthew G. Knepley } 20*6f905325SMatthew G. Knepley ierr = PetscFree(allocated);CHKERRQ(ierr); 21*6f905325SMatthew G. Knepley } 22*6f905325SMatthew G. Knepley ierr = PetscFree(lag->symmetries);CHKERRQ(ierr); 23*6f905325SMatthew G. Knepley } 24*6f905325SMatthew G. Knepley for (i = 0; i < lag->height; i++) { 25*6f905325SMatthew G. Knepley ierr = PetscDualSpaceDestroy(&lag->subspaces[i]);CHKERRQ(ierr); 26*6f905325SMatthew G. Knepley } 27*6f905325SMatthew G. Knepley ierr = PetscFree(lag->subspaces);CHKERRQ(ierr); 28*6f905325SMatthew G. Knepley ierr = PetscFree(lag->numDof);CHKERRQ(ierr); 29*6f905325SMatthew G. Knepley ierr = PetscFree(lag);CHKERRQ(ierr); 30*6f905325SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL);CHKERRQ(ierr); 31*6f905325SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL);CHKERRQ(ierr); 32*6f905325SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", NULL);CHKERRQ(ierr); 33*6f905325SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", NULL);CHKERRQ(ierr); 3420cf1dd8SToby Isaac PetscFunctionReturn(0); 3520cf1dd8SToby Isaac } 3620cf1dd8SToby Isaac 37*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeView_Ascii(PetscDualSpace sp, PetscViewer viewer) 3820cf1dd8SToby Isaac { 3920cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 40*6f905325SMatthew G. Knepley PetscErrorCode ierr; 4120cf1dd8SToby Isaac 4220cf1dd8SToby Isaac PetscFunctionBegin; 43*6f905325SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%s %sLagrange dual space\n", lag->continuous ? "Continuous" : "Discontinuous", lag->tensorSpace ? "tensor " : "");CHKERRQ(ierr); 4420cf1dd8SToby Isaac PetscFunctionReturn(0); 4520cf1dd8SToby Isaac } 4620cf1dd8SToby Isaac 47*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceView_Lagrange(PetscDualSpace sp, PetscViewer viewer) 4820cf1dd8SToby Isaac { 49*6f905325SMatthew G. Knepley PetscBool iascii; 50*6f905325SMatthew G. Knepley PetscErrorCode ierr; 51*6f905325SMatthew G. Knepley 5220cf1dd8SToby Isaac PetscFunctionBegin; 53*6f905325SMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 54*6f905325SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 55*6f905325SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 56*6f905325SMatthew G. Knepley if (iascii) {ierr = PetscDualSpaceLagrangeView_Ascii(sp, viewer);CHKERRQ(ierr);} 5720cf1dd8SToby Isaac PetscFunctionReturn(0); 5820cf1dd8SToby Isaac } 5920cf1dd8SToby Isaac 60*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSetFromOptions_Lagrange(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp) 6120cf1dd8SToby Isaac { 62*6f905325SMatthew G. Knepley PetscBool continuous, tensor, flg; 63*6f905325SMatthew G. Knepley PetscErrorCode ierr; 64*6f905325SMatthew G. Knepley 65*6f905325SMatthew G. Knepley PetscFunctionBegin; 66*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetContinuity(sp, &continuous);CHKERRQ(ierr); 67*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr); 68*6f905325SMatthew G. Knepley ierr = PetscOptionsHead(PetscOptionsObject,"PetscDualSpace Lagrange Options");CHKERRQ(ierr); 69*6f905325SMatthew G. Knepley ierr = PetscOptionsBool("-petscdualspace_lagrange_continuity", "Flag for continuous element", "PetscDualSpaceLagrangeSetContinuity", continuous, &continuous, &flg);CHKERRQ(ierr); 70*6f905325SMatthew G. Knepley if (flg) {ierr = PetscDualSpaceLagrangeSetContinuity(sp, continuous);CHKERRQ(ierr);} 71*6f905325SMatthew G. Knepley ierr = PetscOptionsBool("-petscdualspace_lagrange_tensor", "Flag for tensor dual space", "PetscDualSpaceLagrangeSetContinuity", tensor, &tensor, &flg);CHKERRQ(ierr); 72*6f905325SMatthew G. Knepley if (flg) {ierr = PetscDualSpaceLagrangeSetTensor(sp, tensor);CHKERRQ(ierr);} 73*6f905325SMatthew G. Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 74*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 75*6f905325SMatthew G. Knepley } 76*6f905325SMatthew G. Knepley 77*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceDuplicate_Lagrange(PetscDualSpace sp, PetscDualSpace *spNew) 78*6f905325SMatthew G. Knepley { 79*6f905325SMatthew G. Knepley PetscInt order, Nc; 80*6f905325SMatthew G. Knepley PetscBool cont, tensor; 81*6f905325SMatthew G. Knepley const char *name; 82*6f905325SMatthew G. Knepley PetscErrorCode ierr; 83*6f905325SMatthew G. Knepley 84*6f905325SMatthew G. Knepley PetscFunctionBegin; 85*6f905325SMatthew G. Knepley ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) sp), spNew);CHKERRQ(ierr); 86*6f905325SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) sp, &name);CHKERRQ(ierr); 87*6f905325SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *spNew, name);CHKERRQ(ierr); 88*6f905325SMatthew G. Knepley ierr = PetscDualSpaceSetType(*spNew, PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 89*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetOrder(sp, &order);CHKERRQ(ierr); 90*6f905325SMatthew G. Knepley ierr = PetscDualSpaceSetOrder(*spNew, order);CHKERRQ(ierr); 91*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 92*6f905325SMatthew G. Knepley ierr = PetscDualSpaceSetNumComponents(*spNew, Nc);CHKERRQ(ierr); 93*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetContinuity(sp, &cont);CHKERRQ(ierr); 94*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeSetContinuity(*spNew, cont);CHKERRQ(ierr); 95*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr); 96*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeSetTensor(*spNew, tensor);CHKERRQ(ierr); 97*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 98*6f905325SMatthew G. Knepley } 99*6f905325SMatthew G. Knepley 100*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceGetDimension_SingleCell_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt *dim) 101*6f905325SMatthew G. Knepley { 102*6f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 103*6f905325SMatthew G. Knepley PetscReal D = 1.0; 104*6f905325SMatthew G. Knepley PetscInt n, d; 105*6f905325SMatthew G. Knepley PetscErrorCode ierr; 106*6f905325SMatthew G. Knepley 107*6f905325SMatthew G. Knepley PetscFunctionBegin; 108*6f905325SMatthew G. Knepley *dim = -1; 109*6f905325SMatthew G. Knepley ierr = DMGetDimension(sp->dm, &n);CHKERRQ(ierr); 110*6f905325SMatthew G. Knepley if (!lag->tensorSpace) { 111*6f905325SMatthew G. Knepley for (d = 1; d <= n; ++d) { 112*6f905325SMatthew G. Knepley D *= ((PetscReal) (order+d))/d; 113*6f905325SMatthew G. Knepley } 114*6f905325SMatthew G. Knepley *dim = (PetscInt) (D + 0.5); 115*6f905325SMatthew G. Knepley } else { 116*6f905325SMatthew G. Knepley *dim = 1; 117*6f905325SMatthew G. Knepley for (d = 0; d < n; ++d) *dim *= (order+1); 118*6f905325SMatthew G. Knepley } 119*6f905325SMatthew G. Knepley *dim *= sp->Nc; 120*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 121*6f905325SMatthew G. Knepley } 122*6f905325SMatthew G. Knepley 123*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceCreateHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp) 124*6f905325SMatthew G. Knepley { 125*6f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 126*6f905325SMatthew G. Knepley PetscBool continuous, tensor; 127*6f905325SMatthew G. Knepley PetscInt order; 128*6f905325SMatthew G. Knepley PetscErrorCode ierr; 129*6f905325SMatthew G. Knepley 130*6f905325SMatthew G. Knepley PetscFunctionBegin; 131*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr); 132*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 133*6f905325SMatthew G. Knepley if (height == 0) { 134*6f905325SMatthew G. Knepley ierr = PetscObjectReference((PetscObject)sp);CHKERRQ(ierr); 135*6f905325SMatthew G. Knepley *bdsp = sp; 136*6f905325SMatthew G. Knepley } else if (continuous == PETSC_FALSE || !order) { 137*6f905325SMatthew G. Knepley *bdsp = NULL; 138*6f905325SMatthew G. Knepley } else { 139*6f905325SMatthew G. Knepley DM dm, K; 140*6f905325SMatthew G. Knepley PetscInt dim; 141*6f905325SMatthew G. Knepley 142*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 143*6f905325SMatthew G. Knepley ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 144*6f905325SMatthew G. Knepley if (height > dim || height < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Asked for dual space at height %d for dimension %d reference element\n",height,dim); 145*6f905325SMatthew G. Knepley ierr = PetscDualSpaceDuplicate(sp,bdsp);CHKERRQ(ierr); 146*6f905325SMatthew G. Knepley ierr = PetscDualSpaceCreateReferenceCell(*bdsp, dim-height, lag->simplexCell, &K);CHKERRQ(ierr); 147*6f905325SMatthew G. Knepley ierr = PetscDualSpaceSetDM(*bdsp, K);CHKERRQ(ierr); 148*6f905325SMatthew G. Knepley ierr = DMDestroy(&K);CHKERRQ(ierr); 149*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetTensor(sp,&tensor);CHKERRQ(ierr); 150*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeSetTensor(*bdsp,tensor);CHKERRQ(ierr); 151*6f905325SMatthew G. Knepley ierr = PetscDualSpaceSetUp(*bdsp);CHKERRQ(ierr); 152*6f905325SMatthew G. Knepley } 153*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 154*6f905325SMatthew G. Knepley } 155*6f905325SMatthew G. Knepley 156*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) 157*6f905325SMatthew G. Knepley { 158*6f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 159*6f905325SMatthew G. Knepley DM dm = sp->dm; 160*6f905325SMatthew G. Knepley PetscInt order = sp->order; 161*6f905325SMatthew G. Knepley PetscInt Nc = sp->Nc; 162*6f905325SMatthew G. Knepley MPI_Comm comm; 163*6f905325SMatthew G. Knepley PetscBool continuous; 164*6f905325SMatthew G. Knepley PetscSection csection; 165*6f905325SMatthew G. Knepley Vec coordinates; 166*6f905325SMatthew G. Knepley PetscReal *qpoints, *qweights; 167*6f905325SMatthew G. Knepley PetscInt depth, dim, pdimMax, pStart, pEnd, p, *pStratStart, *pStratEnd, coneSize, d, f = 0, c; 168*6f905325SMatthew G. Knepley PetscBool simplex, tensorSpace; 169*6f905325SMatthew G. Knepley PetscErrorCode ierr; 170*6f905325SMatthew G. Knepley 171*6f905325SMatthew G. Knepley PetscFunctionBegin; 172*6f905325SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) sp, &comm);CHKERRQ(ierr); 173*6f905325SMatthew G. Knepley if (!order) lag->continuous = PETSC_FALSE; 174*6f905325SMatthew G. Knepley continuous = lag->continuous; 175*6f905325SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 176*6f905325SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 177*6f905325SMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 178*6f905325SMatthew G. Knepley ierr = PetscCalloc1(dim+1, &lag->numDof);CHKERRQ(ierr); 179*6f905325SMatthew G. Knepley ierr = PetscMalloc2(depth+1,&pStratStart,depth+1,&pStratEnd);CHKERRQ(ierr); 180*6f905325SMatthew G. Knepley for (d = 0; d <= depth; ++d) {ierr = DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d]);CHKERRQ(ierr);} 181*6f905325SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, pStratStart[depth], &coneSize);CHKERRQ(ierr); 182*6f905325SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &csection);CHKERRQ(ierr); 183*6f905325SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 184*6f905325SMatthew G. Knepley if (depth == 1) { 185*6f905325SMatthew G. Knepley if (coneSize == dim+1) simplex = PETSC_TRUE; 186*6f905325SMatthew G. Knepley else if (coneSize == 1 << dim) simplex = PETSC_FALSE; 187*6f905325SMatthew G. Knepley else SETERRQ(comm, PETSC_ERR_SUP, "Only support simplices and tensor product cells"); 188*6f905325SMatthew G. Knepley } else if (depth == dim) { 189*6f905325SMatthew G. Knepley if (coneSize == dim+1) simplex = PETSC_TRUE; 190*6f905325SMatthew G. Knepley else if (coneSize == 2 * dim) simplex = PETSC_FALSE; 191*6f905325SMatthew G. Knepley else SETERRQ(comm, PETSC_ERR_SUP, "Only support simplices and tensor product cells"); 192*6f905325SMatthew G. Knepley } else SETERRQ(comm, PETSC_ERR_SUP, "Only support cell-vertex meshes or fully interpolated meshes"); 193*6f905325SMatthew G. Knepley lag->simplexCell = simplex; 194*6f905325SMatthew G. Knepley if (dim > 1 && continuous && lag->simplexCell == lag->tensorSpace) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Mismatching simplex/tensor cells and spaces only allowed for discontinuous elements"); 195*6f905325SMatthew G. Knepley tensorSpace = lag->tensorSpace; 196*6f905325SMatthew G. Knepley lag->height = 0; 197*6f905325SMatthew G. Knepley lag->subspaces = NULL; 198*6f905325SMatthew G. Knepley if (continuous && order > 0 && dim > 0) { 19920cf1dd8SToby Isaac PetscInt i; 20020cf1dd8SToby Isaac 201*6f905325SMatthew G. Knepley lag->height = dim; 202*6f905325SMatthew G. Knepley ierr = PetscMalloc1(dim,&lag->subspaces);CHKERRQ(ierr); 203*6f905325SMatthew G. Knepley ierr = PetscDualSpaceCreateHeightSubspace_Lagrange(sp,1,&lag->subspaces[0]);CHKERRQ(ierr); 204*6f905325SMatthew G. Knepley ierr = PetscDualSpaceSetUp(lag->subspaces[0]);CHKERRQ(ierr); 205*6f905325SMatthew G. Knepley for (i = 1; i < dim; i++) { 206*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetHeightSubspace(lag->subspaces[i-1],1,&lag->subspaces[i]);CHKERRQ(ierr); 207*6f905325SMatthew G. Knepley ierr = PetscObjectReference((PetscObject)(lag->subspaces[i]));CHKERRQ(ierr); 208*6f905325SMatthew G. Knepley } 209*6f905325SMatthew G. Knepley } 210*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, &pdimMax);CHKERRQ(ierr); 211*6f905325SMatthew G. Knepley pdimMax *= (pStratEnd[depth] - pStratStart[depth]); 212*6f905325SMatthew G. Knepley ierr = PetscMalloc1(pdimMax, &sp->functional);CHKERRQ(ierr); 213*6f905325SMatthew G. Knepley if (!dim) { 214*6f905325SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 215*6f905325SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 216*6f905325SMatthew G. Knepley ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 217*6f905325SMatthew G. Knepley ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 218*6f905325SMatthew G. Knepley ierr = PetscQuadratureSetData(sp->functional[f], 0, Nc, 1, NULL, qweights);CHKERRQ(ierr); 219*6f905325SMatthew G. Knepley qweights[c] = 1.0; 220*6f905325SMatthew G. Knepley ++f; 221*6f905325SMatthew G. Knepley lag->numDof[0]++; 222*6f905325SMatthew G. Knepley } 22320cf1dd8SToby Isaac } else { 224*6f905325SMatthew G. Knepley PetscSection section; 225*6f905325SMatthew G. Knepley PetscReal *v0, *hv0, *J, *invJ, detJ, hdetJ; 226*6f905325SMatthew G. Knepley PetscInt *tup; 227*6f905325SMatthew G. Knepley 228*6f905325SMatthew G. Knepley ierr = PetscSectionCreate(PETSC_COMM_SELF,§ion);CHKERRQ(ierr); 229*6f905325SMatthew G. Knepley ierr = PetscSectionSetChart(section,pStart,pEnd);CHKERRQ(ierr); 230*6f905325SMatthew G. Knepley ierr = PetscCalloc5(dim+1,&tup,dim,&v0,dim,&hv0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 231*6f905325SMatthew G. Knepley for (p = pStart; p < pEnd; p++) { 232*6f905325SMatthew G. Knepley PetscInt pointDim, d, nFunc = 0; 233*6f905325SMatthew G. Knepley PetscDualSpace hsp; 234*6f905325SMatthew G. Knepley 235*6f905325SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 236*6f905325SMatthew G. Knepley for (d = 0; d < depth; d++) {if (p >= pStratStart[d] && p < pStratEnd[d]) break;} 237*6f905325SMatthew G. Knepley pointDim = (depth == 1 && d == 1) ? dim : d; 238*6f905325SMatthew G. Knepley hsp = ((pointDim < dim) && lag->subspaces) ? lag->subspaces[dim - pointDim - 1] : NULL; 239*6f905325SMatthew G. Knepley if (hsp) { 240*6f905325SMatthew G. Knepley PetscDualSpace_Lag *hlag = (PetscDualSpace_Lag *) hsp->data; 241*6f905325SMatthew G. Knepley DM hdm; 242*6f905325SMatthew G. Knepley 243*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDM(hsp,&hdm);CHKERRQ(ierr); 244*6f905325SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(hdm, 0, NULL, hv0, NULL, NULL, &hdetJ);CHKERRQ(ierr); 245*6f905325SMatthew G. Knepley nFunc = lag->numDof[pointDim] = hlag->numDof[pointDim]; 246*6f905325SMatthew G. Knepley } 247*6f905325SMatthew G. Knepley if (pointDim == dim) { 248*6f905325SMatthew G. Knepley /* Cells, create for self */ 249*6f905325SMatthew G. Knepley PetscInt orderEff = continuous ? (!tensorSpace ? order-1-dim : order-2) : order; 250*6f905325SMatthew G. Knepley PetscReal denom = continuous ? order : (!tensorSpace ? order+1+dim : order+2); 251*6f905325SMatthew G. Knepley PetscReal numer = (!simplex || !tensorSpace) ? 2. : (2./dim); 252*6f905325SMatthew G. Knepley PetscReal dx = numer/denom; 253*6f905325SMatthew G. Knepley PetscInt cdim, d, d2; 254*6f905325SMatthew G. Knepley 255*6f905325SMatthew G. Knepley if (orderEff < 0) continue; 256*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, orderEff, &cdim);CHKERRQ(ierr); 257*6f905325SMatthew G. Knepley ierr = PetscMemzero(tup,(dim+1)*sizeof(PetscInt));CHKERRQ(ierr); 258*6f905325SMatthew G. Knepley if (!tensorSpace) { 259*6f905325SMatthew G. Knepley while (!tup[dim]) { 260*6f905325SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 261*6f905325SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 262*6f905325SMatthew G. Knepley ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr); 263*6f905325SMatthew G. Knepley ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 264*6f905325SMatthew G. Knepley ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 265*6f905325SMatthew G. Knepley ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr); 266*6f905325SMatthew G. Knepley for (d = 0; d < dim; ++d) { 267*6f905325SMatthew G. Knepley qpoints[d] = v0[d]; 268*6f905325SMatthew G. Knepley for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx); 269*6f905325SMatthew G. Knepley } 270*6f905325SMatthew G. Knepley qweights[c] = 1.0; 271*6f905325SMatthew G. Knepley ++f; 272*6f905325SMatthew G. Knepley } 273*6f905325SMatthew G. Knepley ierr = PetscDualSpaceLatticePointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr); 274*6f905325SMatthew G. Knepley } 275*6f905325SMatthew G. Knepley } else { 276*6f905325SMatthew G. Knepley while (!tup[dim]) { 277*6f905325SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 278*6f905325SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 279*6f905325SMatthew G. Knepley ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr); 280*6f905325SMatthew G. Knepley ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 281*6f905325SMatthew G. Knepley ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 282*6f905325SMatthew G. Knepley ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr); 283*6f905325SMatthew G. Knepley for (d = 0; d < dim; ++d) { 284*6f905325SMatthew G. Knepley qpoints[d] = v0[d]; 285*6f905325SMatthew G. Knepley for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx); 286*6f905325SMatthew G. Knepley } 287*6f905325SMatthew G. Knepley qweights[c] = 1.0; 288*6f905325SMatthew G. Knepley ++f; 289*6f905325SMatthew G. Knepley } 290*6f905325SMatthew G. Knepley ierr = PetscDualSpaceTensorPointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr); 29120cf1dd8SToby Isaac } 29220cf1dd8SToby Isaac } 293*6f905325SMatthew G. Knepley lag->numDof[dim] = cdim; 294*6f905325SMatthew G. Knepley } else { /* transform functionals from subspaces */ 295*6f905325SMatthew G. Knepley PetscInt q; 296*6f905325SMatthew G. Knepley 297*6f905325SMatthew G. Knepley for (q = 0; q < nFunc; q++, f++) { 298*6f905325SMatthew G. Knepley PetscQuadrature fn; 299*6f905325SMatthew G. Knepley PetscInt fdim, Nc, c, nPoints, i; 300*6f905325SMatthew G. Knepley const PetscReal *points; 301*6f905325SMatthew G. Knepley const PetscReal *weights; 302*6f905325SMatthew G. Knepley PetscReal *qpoints; 303*6f905325SMatthew G. Knepley PetscReal *qweights; 304*6f905325SMatthew G. Knepley 305*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(hsp, q, &fn);CHKERRQ(ierr); 306*6f905325SMatthew G. Knepley ierr = PetscQuadratureGetData(fn,&fdim,&Nc,&nPoints,&points,&weights);CHKERRQ(ierr); 307*6f905325SMatthew G. Knepley if (fdim != pointDim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected height dual space dim %D, got %D",pointDim,fdim); 308*6f905325SMatthew G. Knepley ierr = PetscMalloc1(nPoints * dim, &qpoints);CHKERRQ(ierr); 309*6f905325SMatthew G. Knepley ierr = PetscCalloc1(nPoints * Nc, &qweights);CHKERRQ(ierr); 310*6f905325SMatthew G. Knepley for (i = 0; i < nPoints; i++) { 311*6f905325SMatthew G. Knepley PetscInt j, k; 312*6f905325SMatthew G. Knepley PetscReal *qp = &qpoints[i * dim]; 313*6f905325SMatthew G. Knepley 314*6f905325SMatthew G. Knepley for (c = 0; c < Nc; ++c) qweights[i*Nc+c] = weights[i*Nc+c]; 315*6f905325SMatthew G. Knepley for (j = 0; j < dim; ++j) qp[j] = v0[j]; 316*6f905325SMatthew G. Knepley for (j = 0; j < dim; ++j) { 317*6f905325SMatthew G. Knepley for (k = 0; k < pointDim; k++) qp[j] += J[dim * j + k] * (points[pointDim * i + k] - hv0[k]); 318*6f905325SMatthew G. Knepley } 319*6f905325SMatthew G. Knepley } 320*6f905325SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 321*6f905325SMatthew G. Knepley ierr = PetscQuadratureSetOrder(sp->functional[f],0);CHKERRQ(ierr); 322*6f905325SMatthew G. Knepley ierr = PetscQuadratureSetData(sp->functional[f],dim,Nc,nPoints,qpoints,qweights);CHKERRQ(ierr); 323*6f905325SMatthew G. Knepley } 324*6f905325SMatthew G. Knepley } 325*6f905325SMatthew G. Knepley ierr = PetscSectionSetDof(section,p,lag->numDof[pointDim]);CHKERRQ(ierr); 326*6f905325SMatthew G. Knepley } 327*6f905325SMatthew G. Knepley ierr = PetscFree5(tup,v0,hv0,J,invJ);CHKERRQ(ierr); 328*6f905325SMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 329*6f905325SMatthew G. Knepley { /* reorder to closure order */ 330*6f905325SMatthew G. Knepley PetscInt *key, count; 331*6f905325SMatthew G. Knepley PetscQuadrature *reorder = NULL; 332*6f905325SMatthew G. Knepley 333*6f905325SMatthew G. Knepley ierr = PetscCalloc1(f,&key);CHKERRQ(ierr); 334*6f905325SMatthew G. Knepley ierr = PetscMalloc1(f*sp->Nc,&reorder);CHKERRQ(ierr); 335*6f905325SMatthew G. Knepley 336*6f905325SMatthew G. Knepley for (p = pStratStart[depth], count = 0; p < pStratEnd[depth]; p++) { 337*6f905325SMatthew G. Knepley PetscInt *closure = NULL, closureSize, c; 338*6f905325SMatthew G. Knepley 339*6f905325SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 340*6f905325SMatthew G. Knepley for (c = 0; c < closureSize; c++) { 341*6f905325SMatthew G. Knepley PetscInt point = closure[2 * c], dof, off, i; 342*6f905325SMatthew G. Knepley 343*6f905325SMatthew G. Knepley ierr = PetscSectionGetDof(section,point,&dof);CHKERRQ(ierr); 344*6f905325SMatthew G. Knepley ierr = PetscSectionGetOffset(section,point,&off);CHKERRQ(ierr); 345*6f905325SMatthew G. Knepley for (i = 0; i < dof; i++) { 346*6f905325SMatthew G. Knepley PetscInt fi = i + off; 347*6f905325SMatthew G. Knepley if (!key[fi]) { 348*6f905325SMatthew G. Knepley key[fi] = 1; 349*6f905325SMatthew G. Knepley reorder[count++] = sp->functional[fi]; 350*6f905325SMatthew G. Knepley } 351*6f905325SMatthew G. Knepley } 352*6f905325SMatthew G. Knepley } 353*6f905325SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 354*6f905325SMatthew G. Knepley } 355*6f905325SMatthew G. Knepley ierr = PetscFree(sp->functional);CHKERRQ(ierr); 356*6f905325SMatthew G. Knepley sp->functional = reorder; 357*6f905325SMatthew G. Knepley ierr = PetscFree(key);CHKERRQ(ierr); 358*6f905325SMatthew G. Knepley } 359*6f905325SMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 360*6f905325SMatthew G. Knepley } 361*6f905325SMatthew G. Knepley if (pStratEnd[depth] == 1 && f != pdimMax) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vectors %D not equal to dimension %D", f, pdimMax); 362*6f905325SMatthew G. Knepley if (f > pdimMax) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vectors %D is greater than max size %D", f, pdimMax); 363*6f905325SMatthew G. Knepley ierr = PetscFree2(pStratStart, pStratEnd);CHKERRQ(ierr); 36420cf1dd8SToby Isaac PetscFunctionReturn(0); 36520cf1dd8SToby Isaac } 36620cf1dd8SToby Isaac 367*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim) 368*6f905325SMatthew G. Knepley { 369*6f905325SMatthew G. Knepley DM K; 370*6f905325SMatthew G. Knepley const PetscInt *numDof; 371*6f905325SMatthew G. Knepley PetscInt spatialDim, Nc, size = 0, d; 372*6f905325SMatthew G. Knepley PetscErrorCode ierr; 373*6f905325SMatthew G. Knepley 374*6f905325SMatthew G. Knepley PetscFunctionBegin; 375*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp, &K);CHKERRQ(ierr); 376*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetNumDof(sp, &numDof);CHKERRQ(ierr); 377*6f905325SMatthew G. Knepley ierr = DMGetDimension(K, &spatialDim);CHKERRQ(ierr); 378*6f905325SMatthew G. Knepley ierr = DMPlexGetHeightStratum(K, 0, NULL, &Nc);CHKERRQ(ierr); 379*6f905325SMatthew G. Knepley if (Nc == 1) {ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, dim);CHKERRQ(ierr); PetscFunctionReturn(0);} 380*6f905325SMatthew G. Knepley for (d = 0; d <= spatialDim; ++d) { 381*6f905325SMatthew G. Knepley PetscInt pStart, pEnd; 382*6f905325SMatthew G. Knepley 383*6f905325SMatthew G. Knepley ierr = DMPlexGetDepthStratum(K, d, &pStart, &pEnd);CHKERRQ(ierr); 384*6f905325SMatthew G. Knepley size += (pEnd-pStart)*numDof[d]; 385*6f905325SMatthew G. Knepley } 386*6f905325SMatthew G. Knepley *dim = size; 387*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 388*6f905325SMatthew G. Knepley } 389*6f905325SMatthew G. Knepley 390*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceGetNumDof_Lagrange(PetscDualSpace sp, const PetscInt **numDof) 391*6f905325SMatthew G. Knepley { 392*6f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 393*6f905325SMatthew G. Knepley 394*6f905325SMatthew G. Knepley PetscFunctionBegin; 395*6f905325SMatthew G. Knepley *numDof = lag->numDof; 396*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 397*6f905325SMatthew G. Knepley } 398*6f905325SMatthew G. Knepley 399*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceGetHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp) 400*6f905325SMatthew G. Knepley { 401*6f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 402*6f905325SMatthew G. Knepley PetscErrorCode ierr; 403*6f905325SMatthew G. Knepley 404*6f905325SMatthew G. Knepley PetscFunctionBegin; 405*6f905325SMatthew G. Knepley if (height == 0) { 406*6f905325SMatthew G. Knepley *bdsp = sp; 407*6f905325SMatthew G. Knepley } else { 408*6f905325SMatthew G. Knepley DM dm; 409*6f905325SMatthew G. Knepley PetscInt dim; 410*6f905325SMatthew G. Knepley 411*6f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 412*6f905325SMatthew G. Knepley ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 413*6f905325SMatthew G. Knepley if (height > dim || height < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Asked for dual space at height %D for dimension %D reference element\n",height,dim); 414*6f905325SMatthew G. Knepley if (height <= lag->height) {*bdsp = lag->subspaces[height-1];} 415*6f905325SMatthew G. Knepley else {*bdsp = NULL;} 416*6f905325SMatthew G. Knepley } 417*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 418*6f905325SMatthew G. Knepley } 41920cf1dd8SToby Isaac 42020cf1dd8SToby Isaac #define BaryIndex(perEdge,a,b,c) (((b)*(2*perEdge+1-(b)))/2)+(c) 42120cf1dd8SToby Isaac 42220cf1dd8SToby Isaac #define CartIndex(perEdge,a,b) (perEdge*(a)+b) 42320cf1dd8SToby Isaac 42420cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceGetSymmetries_Lagrange(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 42520cf1dd8SToby Isaac { 42620cf1dd8SToby Isaac 42720cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 42820cf1dd8SToby Isaac PetscInt dim, order, p, Nc; 42920cf1dd8SToby Isaac PetscErrorCode ierr; 43020cf1dd8SToby Isaac 43120cf1dd8SToby Isaac PetscFunctionBegin; 43220cf1dd8SToby Isaac ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 43320cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumComponents(sp,&Nc);CHKERRQ(ierr); 43420cf1dd8SToby Isaac ierr = DMGetDimension(sp->dm,&dim);CHKERRQ(ierr); 43520cf1dd8SToby Isaac if (!dim || !lag->continuous || order < 3) PetscFunctionReturn(0); 43620cf1dd8SToby Isaac if (dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Lagrange symmetries not implemented for dim = %D > 3",dim); 43720cf1dd8SToby Isaac if (!lag->symmetries) { /* store symmetries */ 43820cf1dd8SToby Isaac PetscDualSpace hsp; 43920cf1dd8SToby Isaac DM K; 44020cf1dd8SToby Isaac PetscInt numPoints = 1, d; 44120cf1dd8SToby Isaac PetscInt numFaces; 44220cf1dd8SToby Isaac PetscInt ***symmetries; 44320cf1dd8SToby Isaac const PetscInt ***hsymmetries; 44420cf1dd8SToby Isaac 44520cf1dd8SToby Isaac if (lag->simplexCell) { 44620cf1dd8SToby Isaac numFaces = 1 + dim; 44720cf1dd8SToby Isaac for (d = 0; d < dim; d++) numPoints = numPoints * 2 + 1; 448*6f905325SMatthew G. Knepley } else { 44920cf1dd8SToby Isaac numPoints = PetscPowInt(3,dim); 45020cf1dd8SToby Isaac numFaces = 2 * dim; 45120cf1dd8SToby Isaac } 45220cf1dd8SToby Isaac ierr = PetscCalloc1(numPoints,&symmetries);CHKERRQ(ierr); 45320cf1dd8SToby Isaac if (0 < dim && dim < 3) { /* compute self symmetries */ 45420cf1dd8SToby Isaac PetscInt **cellSymmetries; 45520cf1dd8SToby Isaac 45620cf1dd8SToby Isaac lag->numSelfSym = 2 * numFaces; 45720cf1dd8SToby Isaac lag->selfSymOff = numFaces; 45820cf1dd8SToby Isaac ierr = PetscCalloc1(2*numFaces,&cellSymmetries);CHKERRQ(ierr); 45920cf1dd8SToby Isaac /* we want to be able to index symmetries directly with the orientations, which range from [-numFaces,numFaces) */ 46020cf1dd8SToby Isaac symmetries[0] = &cellSymmetries[numFaces]; 46120cf1dd8SToby Isaac if (dim == 1) { 46220cf1dd8SToby Isaac PetscInt dofPerEdge = order - 1; 46320cf1dd8SToby Isaac 46420cf1dd8SToby Isaac if (dofPerEdge > 1) { 46520cf1dd8SToby Isaac PetscInt i, j, *reverse; 46620cf1dd8SToby Isaac 46720cf1dd8SToby Isaac ierr = PetscMalloc1(dofPerEdge*Nc,&reverse);CHKERRQ(ierr); 46820cf1dd8SToby Isaac for (i = 0; i < dofPerEdge; i++) { 46920cf1dd8SToby Isaac for (j = 0; j < Nc; j++) { 47020cf1dd8SToby Isaac reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j; 47120cf1dd8SToby Isaac } 47220cf1dd8SToby Isaac } 47320cf1dd8SToby Isaac symmetries[0][-2] = reverse; 47420cf1dd8SToby Isaac 47520cf1dd8SToby Isaac /* yes, this is redundant, but it makes it easier to cleanup if I don't have to worry about what not to free */ 47620cf1dd8SToby Isaac ierr = PetscMalloc1(dofPerEdge*Nc,&reverse);CHKERRQ(ierr); 47720cf1dd8SToby Isaac for (i = 0; i < dofPerEdge; i++) { 47820cf1dd8SToby Isaac for (j = 0; j < Nc; j++) { 47920cf1dd8SToby Isaac reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j; 48020cf1dd8SToby Isaac } 48120cf1dd8SToby Isaac } 48220cf1dd8SToby Isaac symmetries[0][1] = reverse; 48320cf1dd8SToby Isaac } 48420cf1dd8SToby Isaac } else { 48520cf1dd8SToby Isaac PetscInt dofPerEdge = lag->simplexCell ? (order - 2) : (order - 1), s; 48620cf1dd8SToby Isaac PetscInt dofPerFace; 48720cf1dd8SToby Isaac 48820cf1dd8SToby Isaac if (dofPerEdge > 1) { 48920cf1dd8SToby Isaac for (s = -numFaces; s < numFaces; s++) { 49020cf1dd8SToby Isaac PetscInt *sym, i, j, k, l; 49120cf1dd8SToby Isaac 49220cf1dd8SToby Isaac if (!s) continue; 49320cf1dd8SToby Isaac if (lag->simplexCell) { 49420cf1dd8SToby Isaac dofPerFace = (dofPerEdge * (dofPerEdge + 1))/2; 49520cf1dd8SToby Isaac ierr = PetscMalloc1(Nc*dofPerFace,&sym);CHKERRQ(ierr); 49620cf1dd8SToby Isaac for (j = 0, l = 0; j < dofPerEdge; j++) { 49720cf1dd8SToby Isaac for (k = 0; k < dofPerEdge - j; k++, l++) { 49820cf1dd8SToby Isaac i = dofPerEdge - 1 - j - k; 49920cf1dd8SToby Isaac switch (s) { 50020cf1dd8SToby Isaac case -3: 50120cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,i,k,j); 50220cf1dd8SToby Isaac break; 50320cf1dd8SToby Isaac case -2: 50420cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,j,i,k); 50520cf1dd8SToby Isaac break; 50620cf1dd8SToby Isaac case -1: 50720cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,k,j,i); 50820cf1dd8SToby Isaac break; 50920cf1dd8SToby Isaac case 1: 51020cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,k,i,j); 51120cf1dd8SToby Isaac break; 51220cf1dd8SToby Isaac case 2: 51320cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,j,k,i); 51420cf1dd8SToby Isaac break; 51520cf1dd8SToby Isaac } 51620cf1dd8SToby Isaac } 51720cf1dd8SToby Isaac } 51820cf1dd8SToby Isaac } else { 51920cf1dd8SToby Isaac dofPerFace = dofPerEdge * dofPerEdge; 52020cf1dd8SToby Isaac ierr = PetscMalloc1(Nc*dofPerFace,&sym);CHKERRQ(ierr); 52120cf1dd8SToby Isaac for (j = 0, l = 0; j < dofPerEdge; j++) { 52220cf1dd8SToby Isaac for (k = 0; k < dofPerEdge; k++, l++) { 52320cf1dd8SToby Isaac switch (s) { 52420cf1dd8SToby Isaac case -4: 52520cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,k,j); 52620cf1dd8SToby Isaac break; 52720cf1dd8SToby Isaac case -3: 52820cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),k); 52920cf1dd8SToby Isaac break; 53020cf1dd8SToby Isaac case -2: 53120cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),(dofPerEdge - 1 - j)); 53220cf1dd8SToby Isaac break; 53320cf1dd8SToby Isaac case -1: 53420cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,j,(dofPerEdge - 1 - k)); 53520cf1dd8SToby Isaac break; 53620cf1dd8SToby Isaac case 1: 53720cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),j); 53820cf1dd8SToby Isaac break; 53920cf1dd8SToby Isaac case 2: 54020cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),(dofPerEdge - 1 - k)); 54120cf1dd8SToby Isaac break; 54220cf1dd8SToby Isaac case 3: 54320cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,k,(dofPerEdge - 1 - j)); 54420cf1dd8SToby Isaac break; 54520cf1dd8SToby Isaac } 54620cf1dd8SToby Isaac } 54720cf1dd8SToby Isaac } 54820cf1dd8SToby Isaac } 54920cf1dd8SToby Isaac for (i = 0; i < dofPerFace; i++) { 55020cf1dd8SToby Isaac sym[Nc*i] *= Nc; 55120cf1dd8SToby Isaac for (j = 1; j < Nc; j++) { 55220cf1dd8SToby Isaac sym[Nc*i+j] = sym[Nc*i] + j; 55320cf1dd8SToby Isaac } 55420cf1dd8SToby Isaac } 55520cf1dd8SToby Isaac symmetries[0][s] = sym; 55620cf1dd8SToby Isaac } 55720cf1dd8SToby Isaac } 55820cf1dd8SToby Isaac } 55920cf1dd8SToby Isaac } 56020cf1dd8SToby Isaac ierr = PetscDualSpaceGetHeightSubspace(sp,1,&hsp);CHKERRQ(ierr); 56120cf1dd8SToby Isaac ierr = PetscDualSpaceGetSymmetries(hsp,&hsymmetries,NULL);CHKERRQ(ierr); 56220cf1dd8SToby Isaac if (hsymmetries) { 56320cf1dd8SToby Isaac PetscBool *seen; 56420cf1dd8SToby Isaac const PetscInt *cone; 56520cf1dd8SToby Isaac PetscInt KclosureSize, *Kclosure = NULL; 56620cf1dd8SToby Isaac 56720cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp,&K);CHKERRQ(ierr); 56820cf1dd8SToby Isaac ierr = PetscCalloc1(numPoints,&seen);CHKERRQ(ierr); 56920cf1dd8SToby Isaac ierr = DMPlexGetCone(K,0,&cone);CHKERRQ(ierr); 57020cf1dd8SToby Isaac ierr = DMPlexGetTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);CHKERRQ(ierr); 57120cf1dd8SToby Isaac for (p = 0; p < numFaces; p++) { 57220cf1dd8SToby Isaac PetscInt closureSize, *closure = NULL, q; 57320cf1dd8SToby Isaac 57420cf1dd8SToby Isaac ierr = DMPlexGetTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 57520cf1dd8SToby Isaac for (q = 0; q < closureSize; q++) { 57620cf1dd8SToby Isaac PetscInt point = closure[2*q], r; 57720cf1dd8SToby Isaac 57820cf1dd8SToby Isaac if(!seen[point]) { 57920cf1dd8SToby Isaac for (r = 0; r < KclosureSize; r++) { 58020cf1dd8SToby Isaac if (Kclosure[2 * r] == point) break; 58120cf1dd8SToby Isaac } 58220cf1dd8SToby Isaac seen[point] = PETSC_TRUE; 58320cf1dd8SToby Isaac symmetries[r] = (PetscInt **) hsymmetries[q]; 58420cf1dd8SToby Isaac } 58520cf1dd8SToby Isaac } 58620cf1dd8SToby Isaac ierr = DMPlexRestoreTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 58720cf1dd8SToby Isaac } 58820cf1dd8SToby Isaac ierr = DMPlexRestoreTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);CHKERRQ(ierr); 58920cf1dd8SToby Isaac ierr = PetscFree(seen);CHKERRQ(ierr); 59020cf1dd8SToby Isaac } 59120cf1dd8SToby Isaac lag->symmetries = symmetries; 59220cf1dd8SToby Isaac } 59320cf1dd8SToby Isaac if (perms) *perms = (const PetscInt ***) lag->symmetries; 59420cf1dd8SToby Isaac PetscFunctionReturn(0); 59520cf1dd8SToby Isaac } 59620cf1dd8SToby Isaac 59720cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous) 59820cf1dd8SToby Isaac { 59920cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 60020cf1dd8SToby Isaac 60120cf1dd8SToby Isaac PetscFunctionBegin; 60220cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 60320cf1dd8SToby Isaac PetscValidPointer(continuous, 2); 60420cf1dd8SToby Isaac *continuous = lag->continuous; 60520cf1dd8SToby Isaac PetscFunctionReturn(0); 60620cf1dd8SToby Isaac } 60720cf1dd8SToby Isaac 60820cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous) 60920cf1dd8SToby Isaac { 61020cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 61120cf1dd8SToby Isaac 61220cf1dd8SToby Isaac PetscFunctionBegin; 61320cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 61420cf1dd8SToby Isaac lag->continuous = continuous; 61520cf1dd8SToby Isaac PetscFunctionReturn(0); 61620cf1dd8SToby Isaac } 61720cf1dd8SToby Isaac 61820cf1dd8SToby Isaac /*@ 61920cf1dd8SToby Isaac PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity 62020cf1dd8SToby Isaac 62120cf1dd8SToby Isaac Not Collective 62220cf1dd8SToby Isaac 62320cf1dd8SToby Isaac Input Parameter: 62420cf1dd8SToby Isaac . sp - the PetscDualSpace 62520cf1dd8SToby Isaac 62620cf1dd8SToby Isaac Output Parameter: 62720cf1dd8SToby Isaac . continuous - flag for element continuity 62820cf1dd8SToby Isaac 62920cf1dd8SToby Isaac Level: intermediate 63020cf1dd8SToby Isaac 63120cf1dd8SToby Isaac .keywords: PetscDualSpace, Lagrange, continuous, discontinuous 63220cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeSetContinuity() 63320cf1dd8SToby Isaac @*/ 63420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous) 63520cf1dd8SToby Isaac { 63620cf1dd8SToby Isaac PetscErrorCode ierr; 63720cf1dd8SToby Isaac 63820cf1dd8SToby Isaac PetscFunctionBegin; 63920cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 64020cf1dd8SToby Isaac PetscValidPointer(continuous, 2); 64120cf1dd8SToby Isaac ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace,PetscBool*),(sp,continuous));CHKERRQ(ierr); 64220cf1dd8SToby Isaac PetscFunctionReturn(0); 64320cf1dd8SToby Isaac } 64420cf1dd8SToby Isaac 64520cf1dd8SToby Isaac /*@ 64620cf1dd8SToby Isaac PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous 64720cf1dd8SToby Isaac 64820cf1dd8SToby Isaac Logically Collective on PetscDualSpace 64920cf1dd8SToby Isaac 65020cf1dd8SToby Isaac Input Parameters: 65120cf1dd8SToby Isaac + sp - the PetscDualSpace 65220cf1dd8SToby Isaac - continuous - flag for element continuity 65320cf1dd8SToby Isaac 65420cf1dd8SToby Isaac Options Database: 65520cf1dd8SToby Isaac . -petscdualspace_lagrange_continuity <bool> 65620cf1dd8SToby Isaac 65720cf1dd8SToby Isaac Level: intermediate 65820cf1dd8SToby Isaac 65920cf1dd8SToby Isaac .keywords: PetscDualSpace, Lagrange, continuous, discontinuous 66020cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeGetContinuity() 66120cf1dd8SToby Isaac @*/ 66220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous) 66320cf1dd8SToby Isaac { 66420cf1dd8SToby Isaac PetscErrorCode ierr; 66520cf1dd8SToby Isaac 66620cf1dd8SToby Isaac PetscFunctionBegin; 66720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 66820cf1dd8SToby Isaac PetscValidLogicalCollectiveBool(sp, continuous, 2); 66920cf1dd8SToby Isaac ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace,PetscBool),(sp,continuous));CHKERRQ(ierr); 67020cf1dd8SToby Isaac PetscFunctionReturn(0); 67120cf1dd8SToby Isaac } 67220cf1dd8SToby Isaac 673*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Lagrange(PetscDualSpace sp, PetscBool *tensor) 67420cf1dd8SToby Isaac { 67520cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 676*6f905325SMatthew G. Knepley 677*6f905325SMatthew G. Knepley PetscFunctionBegin; 678*6f905325SMatthew G. Knepley *tensor = lag->tensorSpace; 679*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 680*6f905325SMatthew G. Knepley } 681*6f905325SMatthew G. Knepley 682*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Lagrange(PetscDualSpace sp, PetscBool tensor) 683*6f905325SMatthew G. Knepley { 684*6f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 685*6f905325SMatthew G. Knepley 686*6f905325SMatthew G. Knepley PetscFunctionBegin; 687*6f905325SMatthew G. Knepley lag->tensorSpace = tensor; 688*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 689*6f905325SMatthew G. Knepley } 690*6f905325SMatthew G. Knepley 691*6f905325SMatthew G. Knepley /*@ 692*6f905325SMatthew G. Knepley PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space 693*6f905325SMatthew G. Knepley 694*6f905325SMatthew G. Knepley Not collective 695*6f905325SMatthew G. Knepley 696*6f905325SMatthew G. Knepley Input Parameter: 697*6f905325SMatthew G. Knepley . sp - The PetscDualSpace 698*6f905325SMatthew G. Knepley 699*6f905325SMatthew G. Knepley Output Parameter: 700*6f905325SMatthew G. Knepley . tensor - Whether the dual space has tensor layout (vs. simplicial) 701*6f905325SMatthew G. Knepley 702*6f905325SMatthew G. Knepley Level: intermediate 703*6f905325SMatthew G. Knepley 704*6f905325SMatthew G. Knepley .seealso: PetscDualSpaceLagrangeSetTensor(), PetscDualSpaceCreate() 705*6f905325SMatthew G. Knepley @*/ 706*6f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor) 707*6f905325SMatthew G. Knepley { 70820cf1dd8SToby Isaac PetscErrorCode ierr; 70920cf1dd8SToby Isaac 71020cf1dd8SToby Isaac PetscFunctionBegin; 71120cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 712*6f905325SMatthew G. Knepley PetscValidPointer(tensor, 2); 713*6f905325SMatthew G. Knepley ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTensor_C",(PetscDualSpace,PetscBool *),(sp,tensor));CHKERRQ(ierr); 71420cf1dd8SToby Isaac PetscFunctionReturn(0); 71520cf1dd8SToby Isaac } 71620cf1dd8SToby Isaac 717*6f905325SMatthew G. Knepley /*@ 718*6f905325SMatthew G. Knepley PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space 719*6f905325SMatthew G. Knepley 720*6f905325SMatthew G. Knepley Not collective 721*6f905325SMatthew G. Knepley 722*6f905325SMatthew G. Knepley Input Parameters: 723*6f905325SMatthew G. Knepley + sp - The PetscDualSpace 724*6f905325SMatthew G. Knepley - tensor - Whether the dual space has tensor layout (vs. simplicial) 725*6f905325SMatthew G. Knepley 726*6f905325SMatthew G. Knepley Level: intermediate 727*6f905325SMatthew G. Knepley 728*6f905325SMatthew G. Knepley .seealso: PetscDualSpaceLagrangeGetTensor(), PetscDualSpaceCreate() 729*6f905325SMatthew G. Knepley @*/ 730*6f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor) 731*6f905325SMatthew G. Knepley { 732*6f905325SMatthew G. Knepley PetscErrorCode ierr; 733*6f905325SMatthew G. Knepley 734*6f905325SMatthew G. Knepley PetscFunctionBegin; 735*6f905325SMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 736*6f905325SMatthew G. Knepley ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTensor_C",(PetscDualSpace,PetscBool),(sp,tensor));CHKERRQ(ierr); 737*6f905325SMatthew G. Knepley PetscFunctionReturn(0); 738*6f905325SMatthew G. Knepley } 739*6f905325SMatthew G. Knepley 740*6f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) 74120cf1dd8SToby Isaac { 74220cf1dd8SToby Isaac PetscFunctionBegin; 74320cf1dd8SToby Isaac sp->ops->destroy = PetscDualSpaceDestroy_Lagrange; 744*6f905325SMatthew G. Knepley sp->ops->view = PetscDualSpaceView_Lagrange; 745*6f905325SMatthew G. Knepley sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Lagrange; 74620cf1dd8SToby Isaac sp->ops->duplicate = PetscDualSpaceDuplicate_Lagrange; 747*6f905325SMatthew G. Knepley sp->ops->setup = PetscDualSpaceSetUp_Lagrange; 74820cf1dd8SToby Isaac sp->ops->getdimension = PetscDualSpaceGetDimension_Lagrange; 74920cf1dd8SToby Isaac sp->ops->getnumdof = PetscDualSpaceGetNumDof_Lagrange; 75020cf1dd8SToby Isaac sp->ops->getheightsubspace = PetscDualSpaceGetHeightSubspace_Lagrange; 75120cf1dd8SToby Isaac sp->ops->getsymmetries = PetscDualSpaceGetSymmetries_Lagrange; 75220cf1dd8SToby Isaac sp->ops->apply = PetscDualSpaceApplyDefault; 75320cf1dd8SToby Isaac sp->ops->applyall = PetscDualSpaceApplyAllDefault; 75420cf1dd8SToby Isaac sp->ops->createallpoints = PetscDualSpaceCreateAllPointsDefault; 75520cf1dd8SToby Isaac PetscFunctionReturn(0); 75620cf1dd8SToby Isaac } 75720cf1dd8SToby Isaac 75820cf1dd8SToby Isaac /*MC 75920cf1dd8SToby Isaac PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals 76020cf1dd8SToby Isaac 76120cf1dd8SToby Isaac Level: intermediate 76220cf1dd8SToby Isaac 76320cf1dd8SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 76420cf1dd8SToby Isaac M*/ 76520cf1dd8SToby Isaac 76620cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) 76720cf1dd8SToby Isaac { 76820cf1dd8SToby Isaac PetscDualSpace_Lag *lag; 76920cf1dd8SToby Isaac PetscErrorCode ierr; 77020cf1dd8SToby Isaac 77120cf1dd8SToby Isaac PetscFunctionBegin; 77220cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 77320cf1dd8SToby Isaac ierr = PetscNewLog(sp,&lag);CHKERRQ(ierr); 77420cf1dd8SToby Isaac sp->data = lag; 77520cf1dd8SToby Isaac 77620cf1dd8SToby Isaac lag->numDof = NULL; 77720cf1dd8SToby Isaac lag->simplexCell = PETSC_TRUE; 77820cf1dd8SToby Isaac lag->tensorSpace = PETSC_FALSE; 77920cf1dd8SToby Isaac lag->continuous = PETSC_TRUE; 78020cf1dd8SToby Isaac 78120cf1dd8SToby Isaac ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr); 78220cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange);CHKERRQ(ierr); 78320cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange);CHKERRQ(ierr); 78420cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange);CHKERRQ(ierr); 78520cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange);CHKERRQ(ierr); 78620cf1dd8SToby Isaac PetscFunctionReturn(0); 78720cf1dd8SToby Isaac } 78820cf1dd8SToby Isaac 789