120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac #include <petscdmplex.h> 320cf1dd8SToby Isaac 46f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp) 520cf1dd8SToby Isaac { 620cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 76f905325SMatthew G. Knepley PetscInt i; 86f905325SMatthew G. Knepley PetscErrorCode ierr; 920cf1dd8SToby Isaac 1020cf1dd8SToby Isaac PetscFunctionBegin; 116f905325SMatthew G. Knepley if (lag->symmetries) { 126f905325SMatthew G. Knepley PetscInt **selfSyms = lag->symmetries[0]; 136f905325SMatthew G. Knepley 146f905325SMatthew G. Knepley if (selfSyms) { 156f905325SMatthew G. Knepley PetscInt i, **allocated = &selfSyms[-lag->selfSymOff]; 166f905325SMatthew G. Knepley 176f905325SMatthew G. Knepley for (i = 0; i < lag->numSelfSym; i++) { 186f905325SMatthew G. Knepley ierr = PetscFree(allocated[i]);CHKERRQ(ierr); 196f905325SMatthew G. Knepley } 206f905325SMatthew G. Knepley ierr = PetscFree(allocated);CHKERRQ(ierr); 216f905325SMatthew G. Knepley } 226f905325SMatthew G. Knepley ierr = PetscFree(lag->symmetries);CHKERRQ(ierr); 236f905325SMatthew G. Knepley } 246f905325SMatthew G. Knepley for (i = 0; i < lag->height; i++) { 256f905325SMatthew G. Knepley ierr = PetscDualSpaceDestroy(&lag->subspaces[i]);CHKERRQ(ierr); 266f905325SMatthew G. Knepley } 276f905325SMatthew G. Knepley ierr = PetscFree(lag->subspaces);CHKERRQ(ierr); 286f905325SMatthew G. Knepley ierr = PetscFree(lag);CHKERRQ(ierr); 296f905325SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL);CHKERRQ(ierr); 306f905325SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL);CHKERRQ(ierr); 316f905325SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", NULL);CHKERRQ(ierr); 326f905325SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", NULL);CHKERRQ(ierr); 3320cf1dd8SToby Isaac PetscFunctionReturn(0); 3420cf1dd8SToby Isaac } 3520cf1dd8SToby Isaac 366f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeView_Ascii(PetscDualSpace sp, PetscViewer viewer) 3720cf1dd8SToby Isaac { 3820cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 396f905325SMatthew G. Knepley PetscErrorCode ierr; 4020cf1dd8SToby Isaac 4120cf1dd8SToby Isaac PetscFunctionBegin; 426f905325SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%s %sLagrange dual space\n", lag->continuous ? "Continuous" : "Discontinuous", lag->tensorSpace ? "tensor " : "");CHKERRQ(ierr); 4320cf1dd8SToby Isaac PetscFunctionReturn(0); 4420cf1dd8SToby Isaac } 4520cf1dd8SToby Isaac 466f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceView_Lagrange(PetscDualSpace sp, PetscViewer viewer) 4720cf1dd8SToby Isaac { 486f905325SMatthew G. Knepley PetscBool iascii; 496f905325SMatthew G. Knepley PetscErrorCode ierr; 506f905325SMatthew G. Knepley 5120cf1dd8SToby Isaac PetscFunctionBegin; 526f905325SMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 536f905325SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 546f905325SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 556f905325SMatthew G. Knepley if (iascii) {ierr = PetscDualSpaceLagrangeView_Ascii(sp, viewer);CHKERRQ(ierr);} 5620cf1dd8SToby Isaac PetscFunctionReturn(0); 5720cf1dd8SToby Isaac } 5820cf1dd8SToby Isaac 596f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSetFromOptions_Lagrange(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp) 6020cf1dd8SToby Isaac { 616f905325SMatthew G. Knepley PetscBool continuous, tensor, flg; 626f905325SMatthew G. Knepley PetscErrorCode ierr; 636f905325SMatthew G. Knepley 646f905325SMatthew G. Knepley PetscFunctionBegin; 656f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetContinuity(sp, &continuous);CHKERRQ(ierr); 666f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr); 676f905325SMatthew G. Knepley ierr = PetscOptionsHead(PetscOptionsObject,"PetscDualSpace Lagrange Options");CHKERRQ(ierr); 686f905325SMatthew G. Knepley ierr = PetscOptionsBool("-petscdualspace_lagrange_continuity", "Flag for continuous element", "PetscDualSpaceLagrangeSetContinuity", continuous, &continuous, &flg);CHKERRQ(ierr); 696f905325SMatthew G. Knepley if (flg) {ierr = PetscDualSpaceLagrangeSetContinuity(sp, continuous);CHKERRQ(ierr);} 706f905325SMatthew G. Knepley ierr = PetscOptionsBool("-petscdualspace_lagrange_tensor", "Flag for tensor dual space", "PetscDualSpaceLagrangeSetContinuity", tensor, &tensor, &flg);CHKERRQ(ierr); 716f905325SMatthew G. Knepley if (flg) {ierr = PetscDualSpaceLagrangeSetTensor(sp, tensor);CHKERRQ(ierr);} 726f905325SMatthew G. Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 736f905325SMatthew G. Knepley PetscFunctionReturn(0); 746f905325SMatthew G. Knepley } 756f905325SMatthew G. Knepley 76*b4457527SToby Isaac static PetscErrorCode PetscDualSpaceDuplicate_Lagrange(PetscDualSpace sp, PetscDualSpace spNew) 776f905325SMatthew G. Knepley { 786f905325SMatthew G. Knepley PetscBool cont, tensor; 796f905325SMatthew G. Knepley PetscErrorCode ierr; 806f905325SMatthew G. Knepley 816f905325SMatthew G. Knepley PetscFunctionBegin; 826f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetContinuity(sp, &cont);CHKERRQ(ierr); 83*b4457527SToby Isaac ierr = PetscDualSpaceLagrangeSetContinuity(spNew, cont);CHKERRQ(ierr); 846f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr); 85*b4457527SToby Isaac ierr = PetscDualSpaceLagrangeSetTensor(spNew, tensor);CHKERRQ(ierr); 866f905325SMatthew G. Knepley PetscFunctionReturn(0); 876f905325SMatthew G. Knepley } 886f905325SMatthew G. Knepley 896f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceGetDimension_SingleCell_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt *dim) 906f905325SMatthew G. Knepley { 916f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 926f905325SMatthew G. Knepley PetscReal D = 1.0; 936f905325SMatthew G. Knepley PetscInt n, d; 946f905325SMatthew G. Knepley PetscErrorCode ierr; 956f905325SMatthew G. Knepley 966f905325SMatthew G. Knepley PetscFunctionBegin; 976f905325SMatthew G. Knepley *dim = -1; 986f905325SMatthew G. Knepley ierr = DMGetDimension(sp->dm, &n);CHKERRQ(ierr); 996f905325SMatthew G. Knepley if (!lag->tensorSpace) { 1006f905325SMatthew G. Knepley for (d = 1; d <= n; ++d) { 1016f905325SMatthew G. Knepley D *= ((PetscReal) (order+d))/d; 1026f905325SMatthew G. Knepley } 1036f905325SMatthew G. Knepley *dim = (PetscInt) (D + 0.5); 1046f905325SMatthew G. Knepley } else { 1056f905325SMatthew G. Knepley *dim = 1; 1066f905325SMatthew G. Knepley for (d = 0; d < n; ++d) *dim *= (order+1); 1076f905325SMatthew G. Knepley } 1086f905325SMatthew G. Knepley *dim *= sp->Nc; 1096f905325SMatthew G. Knepley PetscFunctionReturn(0); 1106f905325SMatthew G. Knepley } 1116f905325SMatthew G. Knepley 1126f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceCreateHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp) 1136f905325SMatthew G. Knepley { 1146f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 1156f905325SMatthew G. Knepley PetscBool continuous, tensor; 1166f905325SMatthew G. Knepley PetscInt order; 1176f905325SMatthew G. Knepley PetscErrorCode ierr; 1186f905325SMatthew G. Knepley 1196f905325SMatthew G. Knepley PetscFunctionBegin; 1206f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr); 1216f905325SMatthew G. Knepley ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 1226f905325SMatthew G. Knepley if (height == 0) { 1236f905325SMatthew G. Knepley ierr = PetscObjectReference((PetscObject)sp);CHKERRQ(ierr); 1246f905325SMatthew G. Knepley *bdsp = sp; 1256f905325SMatthew G. Knepley } else if (continuous == PETSC_FALSE || !order) { 1266f905325SMatthew G. Knepley *bdsp = NULL; 1276f905325SMatthew G. Knepley } else { 1286f905325SMatthew G. Knepley DM dm, K; 1296f905325SMatthew G. Knepley PetscInt dim; 1306f905325SMatthew G. Knepley 1316f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 1326f905325SMatthew G. Knepley ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 1336f905325SMatthew 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); 1346f905325SMatthew G. Knepley ierr = PetscDualSpaceDuplicate(sp,bdsp);CHKERRQ(ierr); 1356f905325SMatthew G. Knepley ierr = PetscDualSpaceCreateReferenceCell(*bdsp, dim-height, lag->simplexCell, &K);CHKERRQ(ierr); 1366f905325SMatthew G. Knepley ierr = PetscDualSpaceSetDM(*bdsp, K);CHKERRQ(ierr); 1376f905325SMatthew G. Knepley ierr = DMDestroy(&K);CHKERRQ(ierr); 1386f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeGetTensor(sp,&tensor);CHKERRQ(ierr); 1396f905325SMatthew G. Knepley ierr = PetscDualSpaceLagrangeSetTensor(*bdsp,tensor);CHKERRQ(ierr); 1406f905325SMatthew G. Knepley ierr = PetscDualSpaceSetUp(*bdsp);CHKERRQ(ierr); 1416f905325SMatthew G. Knepley } 1426f905325SMatthew G. Knepley PetscFunctionReturn(0); 1436f905325SMatthew G. Knepley } 1446f905325SMatthew G. Knepley 1456f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) 1466f905325SMatthew G. Knepley { 1476f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 1486f905325SMatthew G. Knepley DM dm = sp->dm; 1496f905325SMatthew G. Knepley PetscInt order = sp->order; 1506f905325SMatthew G. Knepley PetscInt Nc = sp->Nc; 1516f905325SMatthew G. Knepley MPI_Comm comm; 1526f905325SMatthew G. Knepley PetscBool continuous; 1536f905325SMatthew G. Knepley PetscSection csection; 1546f905325SMatthew G. Knepley Vec coordinates; 1556f905325SMatthew G. Knepley PetscReal *qpoints, *qweights; 1566f905325SMatthew G. Knepley PetscInt depth, dim, pdimMax, pStart, pEnd, p, *pStratStart, *pStratEnd, coneSize, d, f = 0, c; 1576f905325SMatthew G. Knepley PetscBool simplex, tensorSpace; 1586f905325SMatthew G. Knepley PetscErrorCode ierr; 1596f905325SMatthew G. Knepley 1606f905325SMatthew G. Knepley PetscFunctionBegin; 1616f905325SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) sp, &comm);CHKERRQ(ierr); 1626f905325SMatthew G. Knepley if (!order) lag->continuous = PETSC_FALSE; 1636f905325SMatthew G. Knepley continuous = lag->continuous; 1646f905325SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1656f905325SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1666f905325SMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 1676f905325SMatthew G. Knepley ierr = PetscMalloc2(depth+1,&pStratStart,depth+1,&pStratEnd);CHKERRQ(ierr); 1686f905325SMatthew G. Knepley for (d = 0; d <= depth; ++d) {ierr = DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d]);CHKERRQ(ierr);} 1696f905325SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, pStratStart[depth], &coneSize);CHKERRQ(ierr); 1706f905325SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &csection);CHKERRQ(ierr); 1716f905325SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 1726f905325SMatthew G. Knepley if (depth == 1) { 1736f905325SMatthew G. Knepley if (coneSize == dim+1) simplex = PETSC_TRUE; 1746f905325SMatthew G. Knepley else if (coneSize == 1 << dim) simplex = PETSC_FALSE; 1756f905325SMatthew G. Knepley else SETERRQ(comm, PETSC_ERR_SUP, "Only support simplices and tensor product cells"); 1766f905325SMatthew G. Knepley } else if (depth == dim) { 1776f905325SMatthew G. Knepley if (coneSize == dim+1) simplex = PETSC_TRUE; 1786f905325SMatthew G. Knepley else if (coneSize == 2 * dim) simplex = PETSC_FALSE; 1796f905325SMatthew G. Knepley else SETERRQ(comm, PETSC_ERR_SUP, "Only support simplices and tensor product cells"); 1806f905325SMatthew G. Knepley } else SETERRQ(comm, PETSC_ERR_SUP, "Only support cell-vertex meshes or fully interpolated meshes"); 1816f905325SMatthew G. Knepley lag->simplexCell = simplex; 1826f905325SMatthew 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"); 1836f905325SMatthew G. Knepley tensorSpace = lag->tensorSpace; 1846f905325SMatthew G. Knepley lag->height = 0; 1856f905325SMatthew G. Knepley lag->subspaces = NULL; 1866f905325SMatthew G. Knepley if (continuous && order > 0 && dim > 0) { 18720cf1dd8SToby Isaac PetscInt i; 18820cf1dd8SToby Isaac 1896f905325SMatthew G. Knepley lag->height = dim; 1906f905325SMatthew G. Knepley ierr = PetscMalloc1(dim,&lag->subspaces);CHKERRQ(ierr); 1916f905325SMatthew G. Knepley ierr = PetscDualSpaceCreateHeightSubspace_Lagrange(sp,1,&lag->subspaces[0]);CHKERRQ(ierr); 1926f905325SMatthew G. Knepley ierr = PetscDualSpaceSetUp(lag->subspaces[0]);CHKERRQ(ierr); 1936f905325SMatthew G. Knepley for (i = 1; i < dim; i++) { 1946f905325SMatthew G. Knepley ierr = PetscDualSpaceGetHeightSubspace(lag->subspaces[i-1],1,&lag->subspaces[i]);CHKERRQ(ierr); 1956f905325SMatthew G. Knepley ierr = PetscObjectReference((PetscObject)(lag->subspaces[i]));CHKERRQ(ierr); 1966f905325SMatthew G. Knepley } 1976f905325SMatthew G. Knepley } 1986f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, &pdimMax);CHKERRQ(ierr); 1996f905325SMatthew G. Knepley pdimMax *= (pStratEnd[depth] - pStratStart[depth]); 2006f905325SMatthew G. Knepley ierr = PetscMalloc1(pdimMax, &sp->functional);CHKERRQ(ierr); 2016f905325SMatthew G. Knepley if (!dim) { 2026f905325SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 2036f905325SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 2046f905325SMatthew G. Knepley ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 2056f905325SMatthew G. Knepley ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 2066f905325SMatthew G. Knepley ierr = PetscQuadratureSetData(sp->functional[f], 0, Nc, 1, NULL, qweights);CHKERRQ(ierr); 2076f905325SMatthew G. Knepley qweights[c] = 1.0; 2086f905325SMatthew G. Knepley ++f; 2096f905325SMatthew G. Knepley } 21020cf1dd8SToby Isaac } else { 2116f905325SMatthew G. Knepley PetscSection section; 2126f905325SMatthew G. Knepley PetscReal *v0, *hv0, *J, *invJ, detJ, hdetJ; 2136f905325SMatthew G. Knepley PetscInt *tup; 2146f905325SMatthew G. Knepley 2156f905325SMatthew G. Knepley ierr = PetscSectionCreate(PETSC_COMM_SELF,§ion);CHKERRQ(ierr); 2166f905325SMatthew G. Knepley ierr = PetscSectionSetChart(section,pStart,pEnd);CHKERRQ(ierr); 2176f905325SMatthew G. Knepley ierr = PetscCalloc5(dim+1,&tup,dim,&v0,dim,&hv0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 2186f905325SMatthew G. Knepley for (p = pStart; p < pEnd; p++) { 2196f905325SMatthew G. Knepley PetscInt pointDim, d, nFunc = 0; 2206f905325SMatthew G. Knepley PetscDualSpace hsp; 2216f905325SMatthew G. Knepley 2226f905325SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 2236f905325SMatthew G. Knepley for (d = 0; d < depth; d++) {if (p >= pStratStart[d] && p < pStratEnd[d]) break;} 2246f905325SMatthew G. Knepley pointDim = (depth == 1 && d == 1) ? dim : d; 2256f905325SMatthew G. Knepley hsp = ((pointDim < dim) && lag->subspaces) ? lag->subspaces[dim - pointDim - 1] : NULL; 2266f905325SMatthew G. Knepley if (hsp) { 2276f905325SMatthew G. Knepley DM hdm; 2286f905325SMatthew G. Knepley 2296f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDM(hsp,&hdm);CHKERRQ(ierr); 2306f905325SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(hdm, 0, NULL, hv0, NULL, NULL, &hdetJ);CHKERRQ(ierr); 231*b4457527SToby Isaac nFunc = sp->numDof[pointDim]; 2326f905325SMatthew G. Knepley } 2336f905325SMatthew G. Knepley if (pointDim == dim) { 2346f905325SMatthew G. Knepley /* Cells, create for self */ 2356f905325SMatthew G. Knepley PetscInt orderEff = continuous ? (!tensorSpace ? order-1-dim : order-2) : order; 2366f905325SMatthew G. Knepley PetscReal denom = continuous ? order : (!tensorSpace ? order+1+dim : order+2); 2376f905325SMatthew G. Knepley PetscReal numer = (!simplex || !tensorSpace) ? 2. : (2./dim); 2386f905325SMatthew G. Knepley PetscReal dx = numer/denom; 2396f905325SMatthew G. Knepley PetscInt cdim, d, d2; 2406f905325SMatthew G. Knepley 2416f905325SMatthew G. Knepley if (orderEff < 0) continue; 2426f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, orderEff, &cdim);CHKERRQ(ierr); 243580bdb30SBarry Smith ierr = PetscArrayzero(tup,dim+1);CHKERRQ(ierr); 2446f905325SMatthew G. Knepley if (!tensorSpace) { 2456f905325SMatthew G. Knepley while (!tup[dim]) { 2466f905325SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 2476f905325SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 2486f905325SMatthew G. Knepley ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr); 2496f905325SMatthew G. Knepley ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 2506f905325SMatthew G. Knepley ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 2516f905325SMatthew G. Knepley ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr); 2526f905325SMatthew G. Knepley for (d = 0; d < dim; ++d) { 2536f905325SMatthew G. Knepley qpoints[d] = v0[d]; 2546f905325SMatthew G. Knepley for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx); 2556f905325SMatthew G. Knepley } 2566f905325SMatthew G. Knepley qweights[c] = 1.0; 2576f905325SMatthew G. Knepley ++f; 2586f905325SMatthew G. Knepley } 2596f905325SMatthew G. Knepley ierr = PetscDualSpaceLatticePointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr); 2606f905325SMatthew G. Knepley } 2616f905325SMatthew G. Knepley } else { 2626f905325SMatthew G. Knepley while (!tup[dim]) { 2636f905325SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 2646f905325SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 2656f905325SMatthew G. Knepley ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr); 2666f905325SMatthew G. Knepley ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr); 2676f905325SMatthew G. Knepley ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr); 2686f905325SMatthew G. Knepley ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr); 2696f905325SMatthew G. Knepley for (d = 0; d < dim; ++d) { 2706f905325SMatthew G. Knepley qpoints[d] = v0[d]; 2716f905325SMatthew G. Knepley for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx); 2726f905325SMatthew G. Knepley } 2736f905325SMatthew G. Knepley qweights[c] = 1.0; 2746f905325SMatthew G. Knepley ++f; 2756f905325SMatthew G. Knepley } 2766f905325SMatthew G. Knepley ierr = PetscDualSpaceTensorPointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr); 27720cf1dd8SToby Isaac } 27820cf1dd8SToby Isaac } 2796f905325SMatthew G. Knepley } else { /* transform functionals from subspaces */ 2806f905325SMatthew G. Knepley PetscInt q; 2816f905325SMatthew G. Knepley 2826f905325SMatthew G. Knepley for (q = 0; q < nFunc; q++, f++) { 2836f905325SMatthew G. Knepley PetscQuadrature fn; 2846f905325SMatthew G. Knepley PetscInt fdim, Nc, c, nPoints, i; 2856f905325SMatthew G. Knepley const PetscReal *points; 2866f905325SMatthew G. Knepley const PetscReal *weights; 2876f905325SMatthew G. Knepley PetscReal *qpoints; 2886f905325SMatthew G. Knepley PetscReal *qweights; 2896f905325SMatthew G. Knepley 2906f905325SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(hsp, q, &fn);CHKERRQ(ierr); 2916f905325SMatthew G. Knepley ierr = PetscQuadratureGetData(fn,&fdim,&Nc,&nPoints,&points,&weights);CHKERRQ(ierr); 2926f905325SMatthew G. Knepley if (fdim != pointDim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected height dual space dim %D, got %D",pointDim,fdim); 2936f905325SMatthew G. Knepley ierr = PetscMalloc1(nPoints * dim, &qpoints);CHKERRQ(ierr); 2946f905325SMatthew G. Knepley ierr = PetscCalloc1(nPoints * Nc, &qweights);CHKERRQ(ierr); 2956f905325SMatthew G. Knepley for (i = 0; i < nPoints; i++) { 2966f905325SMatthew G. Knepley PetscInt j, k; 2976f905325SMatthew G. Knepley PetscReal *qp = &qpoints[i * dim]; 2986f905325SMatthew G. Knepley 2996f905325SMatthew G. Knepley for (c = 0; c < Nc; ++c) qweights[i*Nc+c] = weights[i*Nc+c]; 3006f905325SMatthew G. Knepley for (j = 0; j < dim; ++j) qp[j] = v0[j]; 3016f905325SMatthew G. Knepley for (j = 0; j < dim; ++j) { 3026f905325SMatthew G. Knepley for (k = 0; k < pointDim; k++) qp[j] += J[dim * j + k] * (points[pointDim * i + k] - hv0[k]); 3036f905325SMatthew G. Knepley } 3046f905325SMatthew G. Knepley } 3056f905325SMatthew G. Knepley ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr); 3066f905325SMatthew G. Knepley ierr = PetscQuadratureSetOrder(sp->functional[f],0);CHKERRQ(ierr); 3076f905325SMatthew G. Knepley ierr = PetscQuadratureSetData(sp->functional[f],dim,Nc,nPoints,qpoints,qweights);CHKERRQ(ierr); 3086f905325SMatthew G. Knepley } 3096f905325SMatthew G. Knepley } 3106f905325SMatthew G. Knepley } 3116f905325SMatthew G. Knepley ierr = PetscFree5(tup,v0,hv0,J,invJ);CHKERRQ(ierr); 3126f905325SMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 3136f905325SMatthew G. Knepley { /* reorder to closure order */ 3146f905325SMatthew G. Knepley PetscInt *key, count; 3156f905325SMatthew G. Knepley PetscQuadrature *reorder = NULL; 3166f905325SMatthew G. Knepley 3176f905325SMatthew G. Knepley ierr = PetscCalloc1(f,&key);CHKERRQ(ierr); 3186f905325SMatthew G. Knepley ierr = PetscMalloc1(f*sp->Nc,&reorder);CHKERRQ(ierr); 3196f905325SMatthew G. Knepley 3206f905325SMatthew G. Knepley for (p = pStratStart[depth], count = 0; p < pStratEnd[depth]; p++) { 3216f905325SMatthew G. Knepley PetscInt *closure = NULL, closureSize, c; 3226f905325SMatthew G. Knepley 3236f905325SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 3246f905325SMatthew G. Knepley for (c = 0; c < closureSize; c++) { 3256f905325SMatthew G. Knepley PetscInt point = closure[2 * c], dof, off, i; 3266f905325SMatthew G. Knepley 3276f905325SMatthew G. Knepley ierr = PetscSectionGetDof(section,point,&dof);CHKERRQ(ierr); 3286f905325SMatthew G. Knepley ierr = PetscSectionGetOffset(section,point,&off);CHKERRQ(ierr); 3296f905325SMatthew G. Knepley for (i = 0; i < dof; i++) { 3306f905325SMatthew G. Knepley PetscInt fi = i + off; 3316f905325SMatthew G. Knepley if (!key[fi]) { 3326f905325SMatthew G. Knepley key[fi] = 1; 3336f905325SMatthew G. Knepley reorder[count++] = sp->functional[fi]; 3346f905325SMatthew G. Knepley } 3356f905325SMatthew G. Knepley } 3366f905325SMatthew G. Knepley } 3376f905325SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 3386f905325SMatthew G. Knepley } 3396f905325SMatthew G. Knepley ierr = PetscFree(sp->functional);CHKERRQ(ierr); 3406f905325SMatthew G. Knepley sp->functional = reorder; 3416f905325SMatthew G. Knepley ierr = PetscFree(key);CHKERRQ(ierr); 3426f905325SMatthew G. Knepley } 3436f905325SMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 3446f905325SMatthew G. Knepley } 3456f905325SMatthew 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); 3466f905325SMatthew 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); 3476f905325SMatthew G. Knepley ierr = PetscFree2(pStratStart, pStratEnd);CHKERRQ(ierr); 34820cf1dd8SToby Isaac PetscFunctionReturn(0); 34920cf1dd8SToby Isaac } 35020cf1dd8SToby Isaac 3516f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceGetHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp) 3526f905325SMatthew G. Knepley { 3536f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 3546f905325SMatthew G. Knepley PetscErrorCode ierr; 3556f905325SMatthew G. Knepley 3566f905325SMatthew G. Knepley PetscFunctionBegin; 3576f905325SMatthew G. Knepley if (height == 0) { 3586f905325SMatthew G. Knepley *bdsp = sp; 3596f905325SMatthew G. Knepley } else { 3606f905325SMatthew G. Knepley DM dm; 3616f905325SMatthew G. Knepley PetscInt dim; 3626f905325SMatthew G. Knepley 3636f905325SMatthew G. Knepley ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 3646f905325SMatthew G. Knepley ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 3656f905325SMatthew 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); 3666f905325SMatthew G. Knepley if (height <= lag->height) {*bdsp = lag->subspaces[height-1];} 3676f905325SMatthew G. Knepley else {*bdsp = NULL;} 3686f905325SMatthew G. Knepley } 3696f905325SMatthew G. Knepley PetscFunctionReturn(0); 3706f905325SMatthew G. Knepley } 37120cf1dd8SToby Isaac 37220cf1dd8SToby Isaac #define BaryIndex(perEdge,a,b,c) (((b)*(2*perEdge+1-(b)))/2)+(c) 37320cf1dd8SToby Isaac 37420cf1dd8SToby Isaac #define CartIndex(perEdge,a,b) (perEdge*(a)+b) 37520cf1dd8SToby Isaac 37620cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceGetSymmetries_Lagrange(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 37720cf1dd8SToby Isaac { 37820cf1dd8SToby Isaac 37920cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 38020cf1dd8SToby Isaac PetscInt dim, order, p, Nc; 38120cf1dd8SToby Isaac PetscErrorCode ierr; 38220cf1dd8SToby Isaac 38320cf1dd8SToby Isaac PetscFunctionBegin; 38420cf1dd8SToby Isaac ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 38520cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumComponents(sp,&Nc);CHKERRQ(ierr); 38620cf1dd8SToby Isaac ierr = DMGetDimension(sp->dm,&dim);CHKERRQ(ierr); 38720cf1dd8SToby Isaac if (!dim || !lag->continuous || order < 3) PetscFunctionReturn(0); 38820cf1dd8SToby Isaac if (dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Lagrange symmetries not implemented for dim = %D > 3",dim); 38920cf1dd8SToby Isaac if (!lag->symmetries) { /* store symmetries */ 39020cf1dd8SToby Isaac PetscDualSpace hsp; 39120cf1dd8SToby Isaac DM K; 39220cf1dd8SToby Isaac PetscInt numPoints = 1, d; 39320cf1dd8SToby Isaac PetscInt numFaces; 39420cf1dd8SToby Isaac PetscInt ***symmetries; 39520cf1dd8SToby Isaac const PetscInt ***hsymmetries; 39620cf1dd8SToby Isaac 39720cf1dd8SToby Isaac if (lag->simplexCell) { 39820cf1dd8SToby Isaac numFaces = 1 + dim; 39920cf1dd8SToby Isaac for (d = 0; d < dim; d++) numPoints = numPoints * 2 + 1; 4006f905325SMatthew G. Knepley } else { 40120cf1dd8SToby Isaac numPoints = PetscPowInt(3,dim); 40220cf1dd8SToby Isaac numFaces = 2 * dim; 40320cf1dd8SToby Isaac } 40420cf1dd8SToby Isaac ierr = PetscCalloc1(numPoints,&symmetries);CHKERRQ(ierr); 40520cf1dd8SToby Isaac if (0 < dim && dim < 3) { /* compute self symmetries */ 40620cf1dd8SToby Isaac PetscInt **cellSymmetries; 40720cf1dd8SToby Isaac 40820cf1dd8SToby Isaac lag->numSelfSym = 2 * numFaces; 40920cf1dd8SToby Isaac lag->selfSymOff = numFaces; 41020cf1dd8SToby Isaac ierr = PetscCalloc1(2*numFaces,&cellSymmetries);CHKERRQ(ierr); 41120cf1dd8SToby Isaac /* we want to be able to index symmetries directly with the orientations, which range from [-numFaces,numFaces) */ 41220cf1dd8SToby Isaac symmetries[0] = &cellSymmetries[numFaces]; 41320cf1dd8SToby Isaac if (dim == 1) { 41420cf1dd8SToby Isaac PetscInt dofPerEdge = order - 1; 41520cf1dd8SToby Isaac 41620cf1dd8SToby Isaac if (dofPerEdge > 1) { 41720cf1dd8SToby Isaac PetscInt i, j, *reverse; 41820cf1dd8SToby Isaac 41920cf1dd8SToby Isaac ierr = PetscMalloc1(dofPerEdge*Nc,&reverse);CHKERRQ(ierr); 42020cf1dd8SToby Isaac for (i = 0; i < dofPerEdge; i++) { 42120cf1dd8SToby Isaac for (j = 0; j < Nc; j++) { 42220cf1dd8SToby Isaac reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j; 42320cf1dd8SToby Isaac } 42420cf1dd8SToby Isaac } 42520cf1dd8SToby Isaac symmetries[0][-2] = reverse; 42620cf1dd8SToby Isaac 42720cf1dd8SToby Isaac /* yes, this is redundant, but it makes it easier to cleanup if I don't have to worry about what not to free */ 42820cf1dd8SToby Isaac ierr = PetscMalloc1(dofPerEdge*Nc,&reverse);CHKERRQ(ierr); 42920cf1dd8SToby Isaac for (i = 0; i < dofPerEdge; i++) { 43020cf1dd8SToby Isaac for (j = 0; j < Nc; j++) { 43120cf1dd8SToby Isaac reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j; 43220cf1dd8SToby Isaac } 43320cf1dd8SToby Isaac } 43420cf1dd8SToby Isaac symmetries[0][1] = reverse; 43520cf1dd8SToby Isaac } 43620cf1dd8SToby Isaac } else { 43720cf1dd8SToby Isaac PetscInt dofPerEdge = lag->simplexCell ? (order - 2) : (order - 1), s; 43820cf1dd8SToby Isaac PetscInt dofPerFace; 43920cf1dd8SToby Isaac 44020cf1dd8SToby Isaac if (dofPerEdge > 1) { 44120cf1dd8SToby Isaac for (s = -numFaces; s < numFaces; s++) { 44220cf1dd8SToby Isaac PetscInt *sym, i, j, k, l; 44320cf1dd8SToby Isaac 44420cf1dd8SToby Isaac if (!s) continue; 44520cf1dd8SToby Isaac if (lag->simplexCell) { 44620cf1dd8SToby Isaac dofPerFace = (dofPerEdge * (dofPerEdge + 1))/2; 44720cf1dd8SToby Isaac ierr = PetscMalloc1(Nc*dofPerFace,&sym);CHKERRQ(ierr); 44820cf1dd8SToby Isaac for (j = 0, l = 0; j < dofPerEdge; j++) { 44920cf1dd8SToby Isaac for (k = 0; k < dofPerEdge - j; k++, l++) { 45020cf1dd8SToby Isaac i = dofPerEdge - 1 - j - k; 45120cf1dd8SToby Isaac switch (s) { 45220cf1dd8SToby Isaac case -3: 45320cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,i,k,j); 45420cf1dd8SToby Isaac break; 45520cf1dd8SToby Isaac case -2: 45620cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,j,i,k); 45720cf1dd8SToby Isaac break; 45820cf1dd8SToby Isaac case -1: 45920cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,k,j,i); 46020cf1dd8SToby Isaac break; 46120cf1dd8SToby Isaac case 1: 46220cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,k,i,j); 46320cf1dd8SToby Isaac break; 46420cf1dd8SToby Isaac case 2: 46520cf1dd8SToby Isaac sym[Nc*l] = BaryIndex(dofPerEdge,j,k,i); 46620cf1dd8SToby Isaac break; 46720cf1dd8SToby Isaac } 46820cf1dd8SToby Isaac } 46920cf1dd8SToby Isaac } 47020cf1dd8SToby Isaac } else { 47120cf1dd8SToby Isaac dofPerFace = dofPerEdge * dofPerEdge; 47220cf1dd8SToby Isaac ierr = PetscMalloc1(Nc*dofPerFace,&sym);CHKERRQ(ierr); 47320cf1dd8SToby Isaac for (j = 0, l = 0; j < dofPerEdge; j++) { 47420cf1dd8SToby Isaac for (k = 0; k < dofPerEdge; k++, l++) { 47520cf1dd8SToby Isaac switch (s) { 47620cf1dd8SToby Isaac case -4: 47720cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,k,j); 47820cf1dd8SToby Isaac break; 47920cf1dd8SToby Isaac case -3: 48020cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),k); 48120cf1dd8SToby Isaac break; 48220cf1dd8SToby Isaac case -2: 48320cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),(dofPerEdge - 1 - j)); 48420cf1dd8SToby Isaac break; 48520cf1dd8SToby Isaac case -1: 48620cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,j,(dofPerEdge - 1 - k)); 48720cf1dd8SToby Isaac break; 48820cf1dd8SToby Isaac case 1: 48920cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),j); 49020cf1dd8SToby Isaac break; 49120cf1dd8SToby Isaac case 2: 49220cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),(dofPerEdge - 1 - k)); 49320cf1dd8SToby Isaac break; 49420cf1dd8SToby Isaac case 3: 49520cf1dd8SToby Isaac sym[Nc*l] = CartIndex(dofPerEdge,k,(dofPerEdge - 1 - j)); 49620cf1dd8SToby Isaac break; 49720cf1dd8SToby Isaac } 49820cf1dd8SToby Isaac } 49920cf1dd8SToby Isaac } 50020cf1dd8SToby Isaac } 50120cf1dd8SToby Isaac for (i = 0; i < dofPerFace; i++) { 50220cf1dd8SToby Isaac sym[Nc*i] *= Nc; 50320cf1dd8SToby Isaac for (j = 1; j < Nc; j++) { 50420cf1dd8SToby Isaac sym[Nc*i+j] = sym[Nc*i] + j; 50520cf1dd8SToby Isaac } 50620cf1dd8SToby Isaac } 50720cf1dd8SToby Isaac symmetries[0][s] = sym; 50820cf1dd8SToby Isaac } 50920cf1dd8SToby Isaac } 51020cf1dd8SToby Isaac } 51120cf1dd8SToby Isaac } 51220cf1dd8SToby Isaac ierr = PetscDualSpaceGetHeightSubspace(sp,1,&hsp);CHKERRQ(ierr); 51320cf1dd8SToby Isaac ierr = PetscDualSpaceGetSymmetries(hsp,&hsymmetries,NULL);CHKERRQ(ierr); 51420cf1dd8SToby Isaac if (hsymmetries) { 51520cf1dd8SToby Isaac PetscBool *seen; 51620cf1dd8SToby Isaac const PetscInt *cone; 51720cf1dd8SToby Isaac PetscInt KclosureSize, *Kclosure = NULL; 51820cf1dd8SToby Isaac 51920cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp,&K);CHKERRQ(ierr); 52020cf1dd8SToby Isaac ierr = PetscCalloc1(numPoints,&seen);CHKERRQ(ierr); 52120cf1dd8SToby Isaac ierr = DMPlexGetCone(K,0,&cone);CHKERRQ(ierr); 52220cf1dd8SToby Isaac ierr = DMPlexGetTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);CHKERRQ(ierr); 52320cf1dd8SToby Isaac for (p = 0; p < numFaces; p++) { 52420cf1dd8SToby Isaac PetscInt closureSize, *closure = NULL, q; 52520cf1dd8SToby Isaac 52620cf1dd8SToby Isaac ierr = DMPlexGetTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 52720cf1dd8SToby Isaac for (q = 0; q < closureSize; q++) { 52820cf1dd8SToby Isaac PetscInt point = closure[2*q], r; 52920cf1dd8SToby Isaac 53020cf1dd8SToby Isaac if(!seen[point]) { 53120cf1dd8SToby Isaac for (r = 0; r < KclosureSize; r++) { 53220cf1dd8SToby Isaac if (Kclosure[2 * r] == point) break; 53320cf1dd8SToby Isaac } 53420cf1dd8SToby Isaac seen[point] = PETSC_TRUE; 53520cf1dd8SToby Isaac symmetries[r] = (PetscInt **) hsymmetries[q]; 53620cf1dd8SToby Isaac } 53720cf1dd8SToby Isaac } 53820cf1dd8SToby Isaac ierr = DMPlexRestoreTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 53920cf1dd8SToby Isaac } 54020cf1dd8SToby Isaac ierr = DMPlexRestoreTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);CHKERRQ(ierr); 54120cf1dd8SToby Isaac ierr = PetscFree(seen);CHKERRQ(ierr); 54220cf1dd8SToby Isaac } 54320cf1dd8SToby Isaac lag->symmetries = symmetries; 54420cf1dd8SToby Isaac } 54520cf1dd8SToby Isaac if (perms) *perms = (const PetscInt ***) lag->symmetries; 54620cf1dd8SToby Isaac PetscFunctionReturn(0); 54720cf1dd8SToby Isaac } 54820cf1dd8SToby Isaac 54920cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous) 55020cf1dd8SToby Isaac { 55120cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 55220cf1dd8SToby Isaac 55320cf1dd8SToby Isaac PetscFunctionBegin; 55420cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 55520cf1dd8SToby Isaac PetscValidPointer(continuous, 2); 55620cf1dd8SToby Isaac *continuous = lag->continuous; 55720cf1dd8SToby Isaac PetscFunctionReturn(0); 55820cf1dd8SToby Isaac } 55920cf1dd8SToby Isaac 56020cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous) 56120cf1dd8SToby Isaac { 56220cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 56320cf1dd8SToby Isaac 56420cf1dd8SToby Isaac PetscFunctionBegin; 56520cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 56620cf1dd8SToby Isaac lag->continuous = continuous; 56720cf1dd8SToby Isaac PetscFunctionReturn(0); 56820cf1dd8SToby Isaac } 56920cf1dd8SToby Isaac 57020cf1dd8SToby Isaac /*@ 57120cf1dd8SToby Isaac PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity 57220cf1dd8SToby Isaac 57320cf1dd8SToby Isaac Not Collective 57420cf1dd8SToby Isaac 57520cf1dd8SToby Isaac Input Parameter: 57620cf1dd8SToby Isaac . sp - the PetscDualSpace 57720cf1dd8SToby Isaac 57820cf1dd8SToby Isaac Output Parameter: 57920cf1dd8SToby Isaac . continuous - flag for element continuity 58020cf1dd8SToby Isaac 58120cf1dd8SToby Isaac Level: intermediate 58220cf1dd8SToby Isaac 58320cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeSetContinuity() 58420cf1dd8SToby Isaac @*/ 58520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous) 58620cf1dd8SToby Isaac { 58720cf1dd8SToby Isaac PetscErrorCode ierr; 58820cf1dd8SToby Isaac 58920cf1dd8SToby Isaac PetscFunctionBegin; 59020cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 59120cf1dd8SToby Isaac PetscValidPointer(continuous, 2); 59220cf1dd8SToby Isaac ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace,PetscBool*),(sp,continuous));CHKERRQ(ierr); 59320cf1dd8SToby Isaac PetscFunctionReturn(0); 59420cf1dd8SToby Isaac } 59520cf1dd8SToby Isaac 59620cf1dd8SToby Isaac /*@ 59720cf1dd8SToby Isaac PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous 59820cf1dd8SToby Isaac 599d083f849SBarry Smith Logically Collective on sp 60020cf1dd8SToby Isaac 60120cf1dd8SToby Isaac Input Parameters: 60220cf1dd8SToby Isaac + sp - the PetscDualSpace 60320cf1dd8SToby Isaac - continuous - flag for element continuity 60420cf1dd8SToby Isaac 60520cf1dd8SToby Isaac Options Database: 60620cf1dd8SToby Isaac . -petscdualspace_lagrange_continuity <bool> 60720cf1dd8SToby Isaac 60820cf1dd8SToby Isaac Level: intermediate 60920cf1dd8SToby Isaac 61020cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeGetContinuity() 61120cf1dd8SToby Isaac @*/ 61220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous) 61320cf1dd8SToby Isaac { 61420cf1dd8SToby Isaac PetscErrorCode ierr; 61520cf1dd8SToby Isaac 61620cf1dd8SToby Isaac PetscFunctionBegin; 61720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 61820cf1dd8SToby Isaac PetscValidLogicalCollectiveBool(sp, continuous, 2); 61920cf1dd8SToby Isaac ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace,PetscBool),(sp,continuous));CHKERRQ(ierr); 62020cf1dd8SToby Isaac PetscFunctionReturn(0); 62120cf1dd8SToby Isaac } 62220cf1dd8SToby Isaac 6236f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Lagrange(PetscDualSpace sp, PetscBool *tensor) 62420cf1dd8SToby Isaac { 62520cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 6266f905325SMatthew G. Knepley 6276f905325SMatthew G. Knepley PetscFunctionBegin; 6286f905325SMatthew G. Knepley *tensor = lag->tensorSpace; 6296f905325SMatthew G. Knepley PetscFunctionReturn(0); 6306f905325SMatthew G. Knepley } 6316f905325SMatthew G. Knepley 6326f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Lagrange(PetscDualSpace sp, PetscBool tensor) 6336f905325SMatthew G. Knepley { 6346f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 6356f905325SMatthew G. Knepley 6366f905325SMatthew G. Knepley PetscFunctionBegin; 6376f905325SMatthew G. Knepley lag->tensorSpace = tensor; 6386f905325SMatthew G. Knepley PetscFunctionReturn(0); 6396f905325SMatthew G. Knepley } 6406f905325SMatthew G. Knepley 6416f905325SMatthew G. Knepley /*@ 6426f905325SMatthew G. Knepley PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space 6436f905325SMatthew G. Knepley 6446f905325SMatthew G. Knepley Not collective 6456f905325SMatthew G. Knepley 6466f905325SMatthew G. Knepley Input Parameter: 6476f905325SMatthew G. Knepley . sp - The PetscDualSpace 6486f905325SMatthew G. Knepley 6496f905325SMatthew G. Knepley Output Parameter: 6506f905325SMatthew G. Knepley . tensor - Whether the dual space has tensor layout (vs. simplicial) 6516f905325SMatthew G. Knepley 6526f905325SMatthew G. Knepley Level: intermediate 6536f905325SMatthew G. Knepley 6546f905325SMatthew G. Knepley .seealso: PetscDualSpaceLagrangeSetTensor(), PetscDualSpaceCreate() 6556f905325SMatthew G. Knepley @*/ 6566f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor) 6576f905325SMatthew G. Knepley { 65820cf1dd8SToby Isaac PetscErrorCode ierr; 65920cf1dd8SToby Isaac 66020cf1dd8SToby Isaac PetscFunctionBegin; 66120cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 6626f905325SMatthew G. Knepley PetscValidPointer(tensor, 2); 6636f905325SMatthew G. Knepley ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTensor_C",(PetscDualSpace,PetscBool *),(sp,tensor));CHKERRQ(ierr); 66420cf1dd8SToby Isaac PetscFunctionReturn(0); 66520cf1dd8SToby Isaac } 66620cf1dd8SToby Isaac 6676f905325SMatthew G. Knepley /*@ 6686f905325SMatthew G. Knepley PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space 6696f905325SMatthew G. Knepley 6706f905325SMatthew G. Knepley Not collective 6716f905325SMatthew G. Knepley 6726f905325SMatthew G. Knepley Input Parameters: 6736f905325SMatthew G. Knepley + sp - The PetscDualSpace 6746f905325SMatthew G. Knepley - tensor - Whether the dual space has tensor layout (vs. simplicial) 6756f905325SMatthew G. Knepley 6766f905325SMatthew G. Knepley Level: intermediate 6776f905325SMatthew G. Knepley 6786f905325SMatthew G. Knepley .seealso: PetscDualSpaceLagrangeGetTensor(), PetscDualSpaceCreate() 6796f905325SMatthew G. Knepley @*/ 6806f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor) 6816f905325SMatthew G. Knepley { 6826f905325SMatthew G. Knepley PetscErrorCode ierr; 6836f905325SMatthew G. Knepley 6846f905325SMatthew G. Knepley PetscFunctionBegin; 6856f905325SMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 6866f905325SMatthew G. Knepley ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTensor_C",(PetscDualSpace,PetscBool),(sp,tensor));CHKERRQ(ierr); 6876f905325SMatthew G. Knepley PetscFunctionReturn(0); 6886f905325SMatthew G. Knepley } 6896f905325SMatthew G. Knepley 6906f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) 69120cf1dd8SToby Isaac { 69220cf1dd8SToby Isaac PetscFunctionBegin; 69320cf1dd8SToby Isaac sp->ops->destroy = PetscDualSpaceDestroy_Lagrange; 6946f905325SMatthew G. Knepley sp->ops->view = PetscDualSpaceView_Lagrange; 6956f905325SMatthew G. Knepley sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Lagrange; 69620cf1dd8SToby Isaac sp->ops->duplicate = PetscDualSpaceDuplicate_Lagrange; 6976f905325SMatthew G. Knepley sp->ops->setup = PetscDualSpaceSetUp_Lagrange; 698*b4457527SToby Isaac sp->ops->createheightsubspace = PetscDualSpaceGetHeightSubspace_Lagrange; 69920cf1dd8SToby Isaac sp->ops->getsymmetries = PetscDualSpaceGetSymmetries_Lagrange; 70020cf1dd8SToby Isaac sp->ops->apply = PetscDualSpaceApplyDefault; 70120cf1dd8SToby Isaac sp->ops->applyall = PetscDualSpaceApplyAllDefault; 702*b4457527SToby Isaac sp->ops->createalldata = PetscDualSpaceCreateAllDataDefault; 703*b4457527SToby Isaac sp->ops->applyint = PetscDualSpaceApplyInteriorDefault; 704*b4457527SToby Isaac sp->ops->createintdata = PetscDualSpaceCreateInteriorDataDefault; 70520cf1dd8SToby Isaac PetscFunctionReturn(0); 70620cf1dd8SToby Isaac } 70720cf1dd8SToby Isaac 70820cf1dd8SToby Isaac /*MC 70920cf1dd8SToby Isaac PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals 71020cf1dd8SToby Isaac 71120cf1dd8SToby Isaac Level: intermediate 71220cf1dd8SToby Isaac 71320cf1dd8SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 71420cf1dd8SToby Isaac M*/ 71520cf1dd8SToby Isaac 71620cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) 71720cf1dd8SToby Isaac { 71820cf1dd8SToby Isaac PetscDualSpace_Lag *lag; 71920cf1dd8SToby Isaac PetscErrorCode ierr; 72020cf1dd8SToby Isaac 72120cf1dd8SToby Isaac PetscFunctionBegin; 72220cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 72320cf1dd8SToby Isaac ierr = PetscNewLog(sp,&lag);CHKERRQ(ierr); 72420cf1dd8SToby Isaac sp->data = lag; 72520cf1dd8SToby Isaac 72620cf1dd8SToby Isaac lag->simplexCell = PETSC_TRUE; 72720cf1dd8SToby Isaac lag->tensorSpace = PETSC_FALSE; 72820cf1dd8SToby Isaac lag->continuous = PETSC_TRUE; 72920cf1dd8SToby Isaac 73020cf1dd8SToby Isaac ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr); 73120cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange);CHKERRQ(ierr); 73220cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange);CHKERRQ(ierr); 73320cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange);CHKERRQ(ierr); 73420cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange);CHKERRQ(ierr); 73520cf1dd8SToby Isaac PetscFunctionReturn(0); 73620cf1dd8SToby Isaac } 73720cf1dd8SToby Isaac 738