xref: /petsc/src/dm/dt/dualspace/impls/lagrange/tests/ex1.c (revision 5ffaa15ff9d1e1fd99b4639ba3c521328490c7e0)
18a820b8eSToby Isaac 
28a820b8eSToby Isaac #include <petscfe.h>
38a820b8eSToby Isaac #include <petscdmplex.h>
48a820b8eSToby Isaac #include <petsc/private/hashmap.h>
58a820b8eSToby Isaac #include <petsc/private/dmpleximpl.h>
68a820b8eSToby Isaac #include <petsc/private/petscfeimpl.h>
78a820b8eSToby Isaac 
88a820b8eSToby Isaac const char help[] = "Test PETSCDUALSPACELAGRANGE\n";
98a820b8eSToby Isaac 
108a820b8eSToby Isaac typedef struct _PetscHashLagKey
118a820b8eSToby Isaac {
128a820b8eSToby Isaac   PetscInt dim;
138a820b8eSToby Isaac   PetscInt order;
148a820b8eSToby Isaac   PetscInt formDegree;
158a820b8eSToby Isaac   PetscBool trimmed;
168a820b8eSToby Isaac   PetscBool tensor;
178a820b8eSToby Isaac   PetscBool continuous;
188a820b8eSToby Isaac } PetscHashLagKey;
198a820b8eSToby Isaac 
208a820b8eSToby Isaac #define PetscHashLagKeyHash(key) \
218a820b8eSToby Isaac   PetscHashCombine(PetscHashCombine(PetscHashCombine(PetscHashInt((key).dim), \
228a820b8eSToby Isaac                                                      PetscHashInt((key).order)), \
238a820b8eSToby Isaac                                     PetscHashInt((key).formDegree)), \
248a820b8eSToby Isaac                    PetscHashCombine(PetscHashCombine(PetscHashInt((key).trimmed), \
258a820b8eSToby Isaac                                                      PetscHashInt((key).tensor)), \
268a820b8eSToby Isaac                                     PetscHashInt((key).continuous)))
278a820b8eSToby Isaac 
288a820b8eSToby Isaac #define PetscHashLagKeyEqual(k1,k2) \
298a820b8eSToby Isaac   (((k1).dim == (k2).dim) ? \
308a820b8eSToby Isaac    ((k1).order == (k2).order) ? \
318a820b8eSToby Isaac    ((k1).formDegree == (k2).formDegree) ? \
328a820b8eSToby Isaac    ((k1).trimmed == (k2).trimmed) ? \
338a820b8eSToby Isaac    ((k1).tensor == (k2).tensor) ? \
348a820b8eSToby Isaac    ((k1).continuous == (k2).continuous) : 0 : 0 : 0 : 0 : 0)
358a820b8eSToby Isaac 
368a820b8eSToby Isaac PETSC_HASH_MAP(HashLag, PetscHashLagKey, PetscInt, PetscHashLagKeyHash, PetscHashLagKeyEqual, 0)
378a820b8eSToby Isaac 
388a820b8eSToby Isaac static PetscErrorCode ExpectedNumDofs_Total(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs);
398a820b8eSToby Isaac static PetscErrorCode ExpectedNumDofs_Interior(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs);
408a820b8eSToby Isaac 
418a820b8eSToby Isaac static PetscErrorCode ExpectedNumDofs_Total(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs)
428a820b8eSToby Isaac {
438a820b8eSToby Isaac   PetscErrorCode ierr;
448a820b8eSToby Isaac 
458a820b8eSToby Isaac   PetscFunctionBegin;
468a820b8eSToby Isaac   formDegree = PetscAbsInt(formDegree);
478a820b8eSToby Isaac   /* see femtable.org for the source of most of these values */
488a820b8eSToby Isaac   *nDofs = -1;
498a820b8eSToby Isaac   if (tensor == 0) { /* simplex spaces */
508a820b8eSToby Isaac     if (trimmed) {
518a820b8eSToby Isaac       PetscInt rnchooserk;
528a820b8eSToby Isaac       PetscInt rkm1choosek;
538a820b8eSToby Isaac 
548a820b8eSToby Isaac       ierr = PetscDTBinomialInt(order + dim, order + formDegree, &rnchooserk);CHKERRQ(ierr);
558a820b8eSToby Isaac       ierr = PetscDTBinomialInt(order + formDegree - 1, formDegree, &rkm1choosek);CHKERRQ(ierr);
568a820b8eSToby Isaac       *nDofs = rnchooserk * rkm1choosek * nCopies;
578a820b8eSToby Isaac     } else {
588a820b8eSToby Isaac       PetscInt rnchooserk;
598a820b8eSToby Isaac       PetscInt rkchoosek;
608a820b8eSToby Isaac 
618a820b8eSToby Isaac       ierr = PetscDTBinomialInt(order + dim, order + formDegree, &rnchooserk);CHKERRQ(ierr);
628a820b8eSToby Isaac       ierr = PetscDTBinomialInt(order + formDegree, formDegree, &rkchoosek);CHKERRQ(ierr);
638a820b8eSToby Isaac       *nDofs = rnchooserk * rkchoosek * nCopies;
648a820b8eSToby Isaac     }
658a820b8eSToby Isaac   } else if (tensor == 1) { /* hypercubes */
668a820b8eSToby Isaac     if (trimmed) {
678a820b8eSToby Isaac       PetscInt nchoosek;
688a820b8eSToby Isaac       PetscInt rpowk, rp1pownmk;
698a820b8eSToby Isaac 
708a820b8eSToby Isaac       ierr = PetscDTBinomialInt(dim, formDegree, &nchoosek);CHKERRQ(ierr);
718a820b8eSToby Isaac       rpowk = PetscPowInt(order, formDegree);CHKERRQ(ierr);
728a820b8eSToby Isaac       rp1pownmk = PetscPowInt(order + 1, dim - formDegree);CHKERRQ(ierr);
738a820b8eSToby Isaac       *nDofs = nchoosek * rpowk * rp1pownmk * nCopies;
748a820b8eSToby Isaac     } else {
758a820b8eSToby Isaac       PetscInt nchoosek;
768a820b8eSToby Isaac       PetscInt rp1pown;
778a820b8eSToby Isaac 
788a820b8eSToby Isaac       ierr = PetscDTBinomialInt(dim, formDegree, &nchoosek);CHKERRQ(ierr);
798a820b8eSToby Isaac       rp1pown = PetscPowInt(order + 1, dim);CHKERRQ(ierr);
808a820b8eSToby Isaac       *nDofs = nchoosek * rp1pown * nCopies;
818a820b8eSToby Isaac     }
828a820b8eSToby Isaac   } else { /* prism */
838a820b8eSToby Isaac     PetscInt tracek = 0;
848a820b8eSToby Isaac     PetscInt tracekm1 = 0;
858a820b8eSToby Isaac     PetscInt fiber0 = 0;
868a820b8eSToby Isaac     PetscInt fiber1 = 0;
878a820b8eSToby Isaac 
888a820b8eSToby Isaac     if (formDegree < dim) {
898a820b8eSToby Isaac       ierr = ExpectedNumDofs_Total(dim - 1, order, formDegree, trimmed, 0, 1, &tracek);CHKERRQ(ierr);
908a820b8eSToby Isaac       ierr = ExpectedNumDofs_Total(1, order, 0, trimmed, 0, 1, &fiber0);CHKERRQ(ierr);
918a820b8eSToby Isaac     }
928a820b8eSToby Isaac     if (formDegree > 0) {
938a820b8eSToby Isaac       ierr = ExpectedNumDofs_Total(dim - 1, order, formDegree - 1, trimmed, 0, 1, &tracekm1);CHKERRQ(ierr);
948a820b8eSToby Isaac       ierr = ExpectedNumDofs_Total(1, order, 1, trimmed, 0, 1, &fiber1);CHKERRQ(ierr);
958a820b8eSToby Isaac     }
968a820b8eSToby Isaac     *nDofs = (tracek * fiber0 + tracekm1 * fiber1) * nCopies;
978a820b8eSToby Isaac   }
988a820b8eSToby Isaac   PetscFunctionReturn(0);
998a820b8eSToby Isaac }
1008a820b8eSToby Isaac 
1018a820b8eSToby Isaac static PetscErrorCode ExpectedNumDofs_Interior(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed,
1028a820b8eSToby Isaac                                                PetscInt tensor, PetscInt nCopies, PetscInt *nDofs)
1038a820b8eSToby Isaac {
1048a820b8eSToby Isaac   PetscErrorCode ierr;
1058a820b8eSToby Isaac 
1068a820b8eSToby Isaac   PetscFunctionBegin;
1078a820b8eSToby Isaac   formDegree = PetscAbsInt(formDegree);
1088a820b8eSToby Isaac   /* see femtable.org for the source of most of these values */
1098a820b8eSToby Isaac   *nDofs = -1;
1108a820b8eSToby Isaac   if (tensor == 0) { /* simplex spaces */
1118a820b8eSToby Isaac     if (trimmed) {
1128a820b8eSToby Isaac       if (order + formDegree > dim) {
1138a820b8eSToby Isaac         PetscInt eorder = order + formDegree - dim - 1;
1148a820b8eSToby Isaac         PetscInt eformDegree = dim - formDegree;
1158a820b8eSToby Isaac         PetscInt rnchooserk;
1168a820b8eSToby Isaac         PetscInt rkchoosek;
1178a820b8eSToby Isaac 
1188a820b8eSToby Isaac         ierr = PetscDTBinomialInt(eorder + dim, eorder + eformDegree, &rnchooserk);CHKERRQ(ierr);
1198a820b8eSToby Isaac         ierr = PetscDTBinomialInt(eorder + eformDegree, eformDegree, &rkchoosek);CHKERRQ(ierr);
1208a820b8eSToby Isaac         *nDofs = rnchooserk * rkchoosek * nCopies;
1218a820b8eSToby Isaac       } else {
1228a820b8eSToby Isaac         *nDofs = 0;
1238a820b8eSToby Isaac       }
1248a820b8eSToby Isaac 
1258a820b8eSToby Isaac     } else {
1268a820b8eSToby Isaac       if (order + formDegree > dim) {
1278a820b8eSToby Isaac         PetscInt eorder = order + formDegree - dim;
1288a820b8eSToby Isaac         PetscInt eformDegree = dim - formDegree;
1298a820b8eSToby Isaac         PetscInt rnchooserk;
1308a820b8eSToby Isaac         PetscInt rkm1choosek;
1318a820b8eSToby Isaac 
1328a820b8eSToby Isaac         ierr = PetscDTBinomialInt(eorder + dim, eorder + eformDegree, &rnchooserk);CHKERRQ(ierr);
1338a820b8eSToby Isaac         ierr = PetscDTBinomialInt(eorder + eformDegree - 1, eformDegree, &rkm1choosek);CHKERRQ(ierr);
1348a820b8eSToby Isaac         *nDofs = rnchooserk * rkm1choosek * nCopies;
1358a820b8eSToby Isaac       } else {
1368a820b8eSToby Isaac         *nDofs = 0;
1378a820b8eSToby Isaac       }
1388a820b8eSToby Isaac     }
1398a820b8eSToby Isaac   } else if (tensor == 1) { /* hypercubes */
1408a820b8eSToby Isaac     if (dim < 2) {
1418a820b8eSToby Isaac       ierr = ExpectedNumDofs_Interior(dim, order, formDegree, trimmed, 0, nCopies, nDofs);CHKERRQ(ierr);
1428a820b8eSToby Isaac     } else {
1438a820b8eSToby Isaac       PetscInt tracek = 0;
1448a820b8eSToby Isaac       PetscInt tracekm1 = 0;
1458a820b8eSToby Isaac       PetscInt fiber0 = 0;
1468a820b8eSToby Isaac       PetscInt fiber1 = 0;
1478a820b8eSToby Isaac 
1488a820b8eSToby Isaac       if (formDegree < dim) {
1498a820b8eSToby Isaac         ierr = ExpectedNumDofs_Interior(dim - 1, order, formDegree, trimmed, dim > 2 ? 1 : 0, 1, &tracek);CHKERRQ(ierr);
1508a820b8eSToby Isaac         ierr = ExpectedNumDofs_Interior(1, order, 0, trimmed, 0, 1, &fiber0);CHKERRQ(ierr);
1518a820b8eSToby Isaac       }
1528a820b8eSToby Isaac       if (formDegree > 0) {
1538a820b8eSToby Isaac         ierr = ExpectedNumDofs_Interior(dim - 1, order, formDegree - 1, trimmed, dim > 2 ? 1 : 0, 1, &tracekm1);CHKERRQ(ierr);
1548a820b8eSToby Isaac         ierr = ExpectedNumDofs_Interior(1, order, 1, trimmed, 0, 1, &fiber1);CHKERRQ(ierr);
1558a820b8eSToby Isaac       }
1568a820b8eSToby Isaac       *nDofs = (tracek * fiber0 + tracekm1 * fiber1) * nCopies;
1578a820b8eSToby Isaac     }
1588a820b8eSToby Isaac   } else { /* prism */
1598a820b8eSToby Isaac     PetscInt tracek = 0;
1608a820b8eSToby Isaac     PetscInt tracekm1 = 0;
1618a820b8eSToby Isaac     PetscInt fiber0 = 0;
1628a820b8eSToby Isaac     PetscInt fiber1 = 0;
1638a820b8eSToby Isaac 
1648a820b8eSToby Isaac     if (formDegree < dim) {
1658a820b8eSToby Isaac       ierr = ExpectedNumDofs_Interior(dim - 1, order, formDegree, trimmed, 0, 1, &tracek);CHKERRQ(ierr);
1668a820b8eSToby Isaac       ierr = ExpectedNumDofs_Interior(1, order, 0, trimmed, 0, 1, &fiber0);CHKERRQ(ierr);
1678a820b8eSToby Isaac     }
1688a820b8eSToby Isaac     if (formDegree > 0) {
1698a820b8eSToby Isaac       ierr = ExpectedNumDofs_Interior(dim - 1, order, formDegree - 1, trimmed, 0, 1, &tracekm1);CHKERRQ(ierr);
1708a820b8eSToby Isaac       ierr = ExpectedNumDofs_Interior(1, order, 1, trimmed, 0, 1, &fiber1);CHKERRQ(ierr);
1718a820b8eSToby Isaac     }
1728a820b8eSToby Isaac     *nDofs = (tracek * fiber0 + tracekm1 * fiber1) * nCopies;
1738a820b8eSToby Isaac   }
1748a820b8eSToby Isaac   PetscFunctionReturn(0);
1758a820b8eSToby Isaac }
1768a820b8eSToby Isaac 
1778a820b8eSToby Isaac PetscErrorCode testLagrange(PetscHashLag lagTable, DM K, PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscBool continuous, PetscInt nCopies)
1788a820b8eSToby Isaac {
1798a820b8eSToby Isaac   PetscDualSpace  sp;
1808a820b8eSToby Isaac   MPI_Comm        comm = PETSC_COMM_SELF;
1818a820b8eSToby Isaac   PetscInt        Nk;
1828a820b8eSToby Isaac   PetscHashLagKey key;
1838a820b8eSToby Isaac   PetscHashIter   iter;
1848a820b8eSToby Isaac   PetscBool       missing;
1858a820b8eSToby Isaac   PetscInt        spdim, spintdim, exspdim, exspintdim;
1868a820b8eSToby Isaac   PetscErrorCode  ierr;
1878a820b8eSToby Isaac 
1888a820b8eSToby Isaac   PetscFunctionBegin;
1898a820b8eSToby Isaac   ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk);CHKERRQ(ierr);
1908a820b8eSToby Isaac   ierr = PetscDualSpaceCreate(comm, &sp);CHKERRQ(ierr);
1918a820b8eSToby Isaac   ierr = PetscDualSpaceSetType(sp, PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
1928a820b8eSToby Isaac   ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr);
1938a820b8eSToby Isaac   ierr = PetscDualSpaceSetOrder(sp, order);CHKERRQ(ierr);
1948a820b8eSToby Isaac   ierr = PetscDualSpaceSetFormDegree(sp, formDegree);CHKERRQ(ierr);
1958a820b8eSToby Isaac   ierr = PetscDualSpaceSetNumComponents(sp, nCopies * Nk);CHKERRQ(ierr);
1968a820b8eSToby Isaac   ierr = PetscDualSpaceLagrangeSetContinuity(sp, continuous);CHKERRQ(ierr);
1978a820b8eSToby Isaac   ierr = PetscDualSpaceLagrangeSetTensor(sp, (PetscBool) tensor);CHKERRQ(ierr);
1988a820b8eSToby Isaac   ierr = PetscDualSpaceLagrangeSetTrimmed(sp, trimmed);CHKERRQ(ierr);
1998a820b8eSToby Isaac   ierr = PetscInfo7(NULL, "Input: dim %D, order %D, trimmed %D, tensor %D, continuous %D, formDegree %D, nCopies %D\n", dim, order, trimmed, tensor, continuous, formDegree, nCopies);CHKERRQ(ierr);
2008a820b8eSToby Isaac   ierr = ExpectedNumDofs_Total(dim, order, formDegree, trimmed, tensor, nCopies, &exspdim);CHKERRQ(ierr);
2018a820b8eSToby Isaac   if (continuous && dim > 0 && order > 0) {
2028a820b8eSToby Isaac     ierr = ExpectedNumDofs_Interior(dim, order, formDegree, trimmed, tensor, nCopies, &exspintdim);CHKERRQ(ierr);
2038a820b8eSToby Isaac   } else {
2048a820b8eSToby Isaac     exspintdim = exspdim;
2058a820b8eSToby Isaac   }
2068a820b8eSToby Isaac   ierr = PetscDualSpaceSetUp(sp);CHKERRQ(ierr);
2078a820b8eSToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr);
2088a820b8eSToby Isaac   ierr = PetscDualSpaceGetInteriorDimension(sp, &spintdim);CHKERRQ(ierr);
2098a820b8eSToby Isaac   if (spdim != exspdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected dual space dimension %D, got %D\n", exspdim, spdim);CHKERRQ(ierr);
2108a820b8eSToby Isaac   if (spintdim != exspintdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected dual space interior dimension %D, got %D\n", exspintdim, spintdim);CHKERRQ(ierr);
2118a820b8eSToby Isaac   key.dim = dim;
2128a820b8eSToby Isaac   key.formDegree = formDegree;
2138a820b8eSToby Isaac   ierr = PetscDualSpaceGetOrder(sp, &key.order);CHKERRQ(ierr);
2148a820b8eSToby Isaac   ierr = PetscDualSpaceLagrangeGetContinuity(sp, &key.continuous);CHKERRQ(ierr);
2158a820b8eSToby Isaac   if (tensor == 2) {
2168a820b8eSToby Isaac     key.tensor = 2;
2178a820b8eSToby Isaac   } else {
2188a820b8eSToby Isaac     ierr = PetscDualSpaceLagrangeGetTensor(sp, &key.tensor);CHKERRQ(ierr);
2198a820b8eSToby Isaac   }
2208a820b8eSToby Isaac   ierr = PetscDualSpaceLagrangeGetTrimmed(sp, &key.trimmed);CHKERRQ(ierr);
2218a820b8eSToby Isaac   ierr = PetscInfo4(NULL, "After setup:  order %D, trimmed %D, tensor %D, continuous %D\n", key.order, key.trimmed, key.tensor, key.continuous);CHKERRQ(ierr);
2228a820b8eSToby Isaac   ierr = PetscHashLagPut(lagTable, key, &iter, &missing);CHKERRQ(ierr);
2238a820b8eSToby Isaac   if (missing) {
2248a820b8eSToby Isaac     DMPolytopeType type;
2258a820b8eSToby Isaac 
2268a820b8eSToby Isaac     ierr = DMPlexGetCellType(K, 0, &type);CHKERRQ(ierr);
2278a820b8eSToby Isaac     ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "New space: %s, order %D, trimmed %D, tensor %D, continuous %D, form degree %D\n", DMPolytopeTypes[type], order, trimmed, tensor, continuous, formDegree);CHKERRQ(ierr);
2288a820b8eSToby Isaac     ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2298a820b8eSToby Isaac     {
2308a820b8eSToby Isaac       PetscQuadrature intNodes, allNodes;
2318a820b8eSToby Isaac       Mat intMat, allMat;
2328a820b8eSToby Isaac       MatInfo info;
2338a820b8eSToby Isaac       PetscInt i, j, nodeIdxDim, nodeVecDim, nNodes;
2348a820b8eSToby Isaac       const PetscInt *nodeIdx;
2358a820b8eSToby Isaac       const PetscReal *nodeVec;
2368a820b8eSToby Isaac 
2378a820b8eSToby Isaac       PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2388a820b8eSToby Isaac 
2398a820b8eSToby Isaac       ierr = PetscLagNodeIndicesGetData_Internal(lag->allNodeIndices, &nodeIdxDim, &nodeVecDim, &nNodes, &nodeIdx, &nodeVec);CHKERRQ(ierr);
2408a820b8eSToby Isaac       if (nodeVecDim != Nk) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect nodeVecDim");CHKERRQ(ierr);
2418a820b8eSToby Isaac       if (nNodes != spdim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect nNodes");CHKERRQ(ierr);
2428a820b8eSToby Isaac 
2438a820b8eSToby Isaac       ierr = PetscDualSpaceGetAllData(sp, &allNodes, &allMat);CHKERRQ(ierr);
2448a820b8eSToby Isaac 
2458a820b8eSToby Isaac       ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "All nodes:\n");CHKERRQ(ierr);
2468a820b8eSToby Isaac       ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2478a820b8eSToby Isaac       ierr = PetscQuadratureView(allNodes, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2488a820b8eSToby Isaac       ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2498a820b8eSToby Isaac       ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "All node indices:\n");CHKERRQ(ierr);
2508a820b8eSToby Isaac       for (i = 0; i < spdim; i++) {
2518a820b8eSToby Isaac         ierr = PetscPrintf(PETSC_COMM_SELF, "(");CHKERRQ(ierr);
2528a820b8eSToby Isaac         for (j = 0; j < nodeIdxDim; j++) {
2538a820b8eSToby Isaac           ierr = PetscPrintf(PETSC_COMM_SELF, " %D,", nodeIdx[i * nodeIdxDim + j]);CHKERRQ(ierr);
2548a820b8eSToby Isaac         }
2558a820b8eSToby Isaac         ierr = PetscPrintf(PETSC_COMM_SELF, "): [");CHKERRQ(ierr);
2568a820b8eSToby Isaac         for (j = 0; j < nodeVecDim; j++) {
2578a820b8eSToby Isaac           ierr = PetscPrintf(PETSC_COMM_SELF, " %g,", nodeVec[i * nodeVecDim + j]);CHKERRQ(ierr);
2588a820b8eSToby Isaac         }
2598a820b8eSToby Isaac         ierr = PetscPrintf(PETSC_COMM_SELF, "]\n");CHKERRQ(ierr);
2608a820b8eSToby Isaac       }
2618a820b8eSToby Isaac 
2628a820b8eSToby Isaac       ierr = MatGetInfo(allMat, MAT_LOCAL, &info);CHKERRQ(ierr);
2638a820b8eSToby Isaac       ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "All matrix: %D nonzeros\n", (PetscInt) info.nz_used);CHKERRQ(ierr);
2648a820b8eSToby Isaac 
2658a820b8eSToby Isaac       ierr = PetscDualSpaceGetInteriorData(sp, &intNodes, &intMat);CHKERRQ(ierr);
2668a820b8eSToby Isaac       if (intMat && intMat != allMat) {
2678a820b8eSToby Isaac         PetscInt intNodeIdxDim, intNodeVecDim, intNnodes;
2688a820b8eSToby Isaac         const PetscInt *intNodeIdx;
2698a820b8eSToby Isaac         const PetscReal *intNodeVec;
2708a820b8eSToby Isaac         PetscBool same;
2718a820b8eSToby Isaac 
2728a820b8eSToby Isaac         ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior nodes:\n");CHKERRQ(ierr);
2738a820b8eSToby Isaac         ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2748a820b8eSToby Isaac         ierr = PetscQuadratureView(intNodes, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2758a820b8eSToby Isaac         ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
2768a820b8eSToby Isaac 
2778a820b8eSToby Isaac         ierr = MatGetInfo(intMat, MAT_LOCAL, &info);CHKERRQ(ierr);
2788a820b8eSToby Isaac         ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior matrix: %D nonzeros\n", (PetscInt) info.nz_used);CHKERRQ(ierr);
2798a820b8eSToby Isaac         ierr = PetscLagNodeIndicesGetData_Internal(lag->intNodeIndices, &intNodeIdxDim, &intNodeVecDim, &intNnodes, &intNodeIdx, &intNodeVec);CHKERRQ(ierr);
2808a820b8eSToby Isaac         if (intNodeIdxDim != nodeIdxDim || intNodeVecDim != nodeVecDim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node indices not the same shale as all node indices");
2818a820b8eSToby Isaac         if (intNnodes != spintdim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect interior nNodes");CHKERRQ(ierr);
2828a820b8eSToby Isaac         ierr = PetscArraycmp(intNodeIdx, nodeIdx, nodeIdxDim * intNnodes, &same);CHKERRQ(ierr);
2838a820b8eSToby Isaac         if (!same) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node indices not the same as start of all node indices");
2848a820b8eSToby Isaac         ierr = PetscArraycmp(intNodeVec, nodeVec, nodeVecDim * intNnodes, &same);CHKERRQ(ierr);
2858a820b8eSToby Isaac         if (!same) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node vectors not the same as start of all node vectors");
2868a820b8eSToby Isaac       } else if (intMat) {
2878a820b8eSToby Isaac         ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior data is the same as all data\n");CHKERRQ(ierr);
2888a820b8eSToby Isaac         if (intNodes != allNodes) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior nodes should be the same as all nodes");
2898a820b8eSToby Isaac         if (lag->intNodeIndices != lag->allNodeIndices) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node indices should be the same as all node indices");
2908a820b8eSToby Isaac       }
2918a820b8eSToby Isaac     }
2928a820b8eSToby Isaac     if (dim <= 2 && spintdim) {
2938a820b8eSToby Isaac       PetscInt coneSize, o;
2948a820b8eSToby Isaac 
2958a820b8eSToby Isaac       ierr = DMPlexGetConeSize(K, 0, &coneSize);CHKERRQ(ierr);
2968a820b8eSToby Isaac       for (o = -coneSize; o < coneSize; o++) {
2978a820b8eSToby Isaac         Mat symMat;
2988a820b8eSToby Isaac 
2998a820b8eSToby Isaac         ierr = PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(sp, o, &symMat);CHKERRQ(ierr);
3008a820b8eSToby Isaac         ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior node symmetry matrix for orientation %D:\n", o);CHKERRQ(ierr);
3018a820b8eSToby Isaac         ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
3028a820b8eSToby Isaac         ierr = MatView(symMat, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
3038a820b8eSToby Isaac         ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
3048a820b8eSToby Isaac         ierr = MatDestroy(&symMat);CHKERRQ(ierr);
3058a820b8eSToby Isaac       }
3068a820b8eSToby Isaac     }
3078a820b8eSToby Isaac     ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
3088a820b8eSToby Isaac   }
3098a820b8eSToby Isaac   ierr = PetscDualSpaceDestroy(&sp);CHKERRQ(ierr);
3108a820b8eSToby Isaac   PetscFunctionReturn(0);
3118a820b8eSToby Isaac }
3128a820b8eSToby Isaac 
3138a820b8eSToby Isaac static PetscErrorCode DMPlexCreateReferenceWedge(MPI_Comm comm, DM *refdm)
3148a820b8eSToby Isaac {
3158a820b8eSToby Isaac   PetscInt       dim = 3;
3168a820b8eSToby Isaac   DM             rdm;
3178a820b8eSToby Isaac   PetscErrorCode ierr;
3188a820b8eSToby Isaac 
3198a820b8eSToby Isaac   PetscFunctionBeginUser;
3208a820b8eSToby Isaac   ierr = DMCreate(comm, &rdm);CHKERRQ(ierr);
3218a820b8eSToby Isaac   ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr);
3228a820b8eSToby Isaac   ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr);
3238a820b8eSToby Isaac   {
3248a820b8eSToby Isaac     PetscInt    numPoints[4]         = {6, 9, 5, 1};
3258a820b8eSToby Isaac     PetscInt    coneSize[21]         = {5,
3268a820b8eSToby Isaac                                         3, 3,
3278a820b8eSToby Isaac                                         4, 4, 4,
3288a820b8eSToby Isaac                                         2, 2, 2, 2, 2, 2, 2, 2, 2,
3298a820b8eSToby Isaac                                         0, 0, 0, 0, 0, 0};
3308a820b8eSToby Isaac     PetscInt    cones[41]            = {1, 2, 3, 4, 5,
3318a820b8eSToby Isaac                                         6, 7, 8,
3328a820b8eSToby Isaac                                         9, 10, 11,
3338a820b8eSToby Isaac                                         8, 12, 9, 13,
3348a820b8eSToby Isaac                                         7, 14, 10, 12,
3358a820b8eSToby Isaac                                         6, 13, 11, 14,
3368a820b8eSToby Isaac                                         15, 16,  16, 17,  17, 15,
3378a820b8eSToby Isaac                                         18, 19,  19, 20,  20, 18,
3388a820b8eSToby Isaac                                         17, 19,  18, 15,  16, 20};
3398a820b8eSToby Isaac     PetscInt    coneOrientations[41] = {0, 0, 0, 0, 0,
3408a820b8eSToby Isaac                                         0, 0, 0,
3418a820b8eSToby Isaac                                         0, 0, 0,
3428a820b8eSToby Isaac                                         -2,  0, -2,  0,
3438a820b8eSToby Isaac                                         -2,  0, -2, -2,
3448a820b8eSToby Isaac                                         -2, -2, -2, -2,
3458a820b8eSToby Isaac                                         0, 0, 0, 0, 0, 0,
3468a820b8eSToby Isaac                                         0, 0, 0, 0, 0, 0,
3478a820b8eSToby Isaac                                         0, 0, 0, 0, 0, 0};
3488a820b8eSToby Isaac     PetscScalar vertexCoords[18]    = {-1.0, -1.0, -1.0,
3498a820b8eSToby Isaac                                        -1.0,  1.0, -1.0,
3508a820b8eSToby Isaac                                         1.0, -1.0, -1.0,
3518a820b8eSToby Isaac                                        -1.0, -1.0,  1.0,
3528a820b8eSToby Isaac                                         1.0, -1.0,  1.0,
3538a820b8eSToby Isaac                                        -1.0,  1.0,  1.0};
3548a820b8eSToby Isaac 
3558a820b8eSToby Isaac     ierr = DMPlexCreateFromDAG(rdm, 3, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
3568a820b8eSToby Isaac   }
3578a820b8eSToby Isaac   *refdm = rdm;
3588a820b8eSToby Isaac   PetscFunctionReturn(0);
3598a820b8eSToby Isaac }
3608a820b8eSToby Isaac 
3618a820b8eSToby Isaac int main (int argc, char **argv)
3628a820b8eSToby Isaac {
3638a820b8eSToby Isaac   PetscInt        dim;
3648a820b8eSToby Isaac   PetscHashLag    lagTable;
3658a820b8eSToby Isaac   PetscInt        tensorCell;
3668a820b8eSToby Isaac   PetscInt        order, ordermin, ordermax;
3678a820b8eSToby Isaac   PetscBool       continuous;
3688a820b8eSToby Isaac   PetscBool       trimmed;
3698a820b8eSToby Isaac   DM              dm;
3708a820b8eSToby Isaac   PetscErrorCode  ierr;
3718a820b8eSToby Isaac 
3728a820b8eSToby Isaac   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
3738a820b8eSToby Isaac   dim = 3;
3748a820b8eSToby Isaac   tensorCell = 0;
3758a820b8eSToby Isaac   continuous = PETSC_FALSE;
3768a820b8eSToby Isaac   trimmed = PETSC_FALSE;
3778a820b8eSToby Isaac   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for PETSCDUALSPACELAGRANGE test","none");CHKERRQ(ierr);
3788a820b8eSToby Isaac   ierr = PetscOptionsRangeInt("-dim", "The spatial dimension","ex1.c",dim,&dim,NULL,0,3);CHKERRQ(ierr);
379*5ffaa15fSToby Isaac   ierr = PetscOptionsRangeInt("-tensor", "(0) simplex (1) hypercube (2) wedge","ex1.c",tensorCell,&tensorCell,NULL,0,2);CHKERRQ(ierr);
3808a820b8eSToby Isaac   ierr = PetscOptionsBool("-continuous", "Whether the dual space has continuity","ex1.c",continuous,&continuous,NULL);CHKERRQ(ierr);
3818a820b8eSToby Isaac   ierr = PetscOptionsBool("-trimmed", "Whether the dual space matches a trimmed polynomial space","ex1.c",trimmed,&trimmed,NULL);CHKERRQ(ierr);
3828a820b8eSToby Isaac   ierr = PetscOptionsEnd();
3838a820b8eSToby Isaac   ierr = PetscHashLagCreate(&lagTable);CHKERRQ(ierr);
3848a820b8eSToby Isaac 
3858a820b8eSToby Isaac   if (tensorCell < 2) {
3868a820b8eSToby Isaac     ierr = DMPlexCreateReferenceCell(PETSC_COMM_SELF, dim, (PetscBool) !tensorCell, &dm);CHKERRQ(ierr);
3878a820b8eSToby Isaac   } else {
3888a820b8eSToby Isaac     ierr = DMPlexCreateReferenceWedge(PETSC_COMM_SELF, &dm);CHKERRQ(ierr);
3898a820b8eSToby Isaac   }
3908a820b8eSToby Isaac   ordermin = trimmed ? 1 : 0;
3918a820b8eSToby Isaac   ordermax = tensorCell == 2 ? 4 : tensorCell == 1 ? 3 : dim + 2;
3928a820b8eSToby Isaac   for (order = ordermin; order <= ordermax; order++) {
3938a820b8eSToby Isaac     PetscInt formDegree;
3948a820b8eSToby Isaac 
3958a820b8eSToby Isaac     for (formDegree = PetscMin(0,-dim+1); formDegree <= dim; formDegree++) {
3968a820b8eSToby Isaac       PetscInt nCopies;
3978a820b8eSToby Isaac 
3988a820b8eSToby Isaac       for (nCopies = 1; nCopies <= 3; nCopies++) {
3998a820b8eSToby Isaac         ierr = testLagrange(lagTable, dm, dim, order, formDegree, trimmed, (PetscBool) tensorCell, continuous, nCopies);CHKERRQ(ierr);
4008a820b8eSToby Isaac       }
4018a820b8eSToby Isaac     }
4028a820b8eSToby Isaac   }
4038a820b8eSToby Isaac   ierr = DMDestroy(&dm);CHKERRQ(ierr);
4048a820b8eSToby Isaac   ierr = PetscHashLagDestroy(&lagTable);CHKERRQ(ierr);
4058a820b8eSToby Isaac   ierr = PetscFinalize();
4068a820b8eSToby Isaac   return ierr;
4078a820b8eSToby Isaac }
4088a820b8eSToby Isaac 
4098a820b8eSToby Isaac /*TEST
4108a820b8eSToby Isaac 
4118a820b8eSToby Isaac  test:
4128a820b8eSToby Isaac    suffix: 0
4138a820b8eSToby Isaac    args: -dim 0
4148a820b8eSToby Isaac 
4158a820b8eSToby Isaac  test:
4168a820b8eSToby Isaac    suffix: 1_discontinuous_full
4178a820b8eSToby Isaac    args: -dim 1 -continuous 0 -trimmed 0
4188a820b8eSToby Isaac 
4198a820b8eSToby Isaac  test:
4208a820b8eSToby Isaac    suffix: 1_continuous_full
4218a820b8eSToby Isaac    args: -dim 1 -continuous 1 -trimmed 0
4228a820b8eSToby Isaac 
4238a820b8eSToby Isaac  test:
4248a820b8eSToby Isaac    suffix: 2_simplex_discontinuous_full
4258a820b8eSToby Isaac    args: -dim 2 -tensor 0 -continuous 0 -trimmed 0
4268a820b8eSToby Isaac 
4278a820b8eSToby Isaac  test:
4288a820b8eSToby Isaac    suffix: 2_simplex_continuous_full
4298a820b8eSToby Isaac    args: -dim 2 -tensor 0 -continuous 1 -trimmed 0
4308a820b8eSToby Isaac 
4318a820b8eSToby Isaac  test:
4328a820b8eSToby Isaac    suffix: 2_tensor_discontinuous_full
4338a820b8eSToby Isaac    args: -dim 2 -tensor 1 -continuous 0 -trimmed 0
4348a820b8eSToby Isaac 
4358a820b8eSToby Isaac  test:
4368a820b8eSToby Isaac    suffix: 2_tensor_continuous_full
4378a820b8eSToby Isaac    args: -dim 2 -tensor 1 -continuous 1 -trimmed 0
4388a820b8eSToby Isaac 
4398a820b8eSToby Isaac  test:
4408a820b8eSToby Isaac    suffix: 3_simplex_discontinuous_full
4418a820b8eSToby Isaac    args: -dim 3 -tensor 0 -continuous 0 -trimmed 0
4428a820b8eSToby Isaac 
4438a820b8eSToby Isaac  test:
4448a820b8eSToby Isaac    suffix: 3_simplex_continuous_full
4458a820b8eSToby Isaac    args: -dim 3 -tensor 0 -continuous 1 -trimmed 0
4468a820b8eSToby Isaac 
4478a820b8eSToby Isaac  test:
4488a820b8eSToby Isaac    suffix: 3_tensor_discontinuous_full
4498a820b8eSToby Isaac    args: -dim 3 -tensor 1 -continuous 0 -trimmed 0
4508a820b8eSToby Isaac 
4518a820b8eSToby Isaac  test:
4528a820b8eSToby Isaac    suffix: 3_tensor_continuous_full
4538a820b8eSToby Isaac    args: -dim 3 -tensor 1 -continuous 1 -trimmed 0
4548a820b8eSToby Isaac 
4558a820b8eSToby Isaac  test:
4568a820b8eSToby Isaac    suffix: 3_wedge_discontinuous_full
4578a820b8eSToby Isaac    args: -dim 3 -tensor 2 -continuous 0 -trimmed 0
4588a820b8eSToby Isaac 
4598a820b8eSToby Isaac  test:
4608a820b8eSToby Isaac    suffix: 3_wedge_continuous_full
4618a820b8eSToby Isaac    args: -dim 3 -tensor 2 -continuous 1 -trimmed 0
4628a820b8eSToby Isaac 
4638a820b8eSToby Isaac  test:
4648a820b8eSToby Isaac    suffix: 1_discontinuous_trimmed
4658a820b8eSToby Isaac    args: -dim 1 -continuous 0 -trimmed 1
4668a820b8eSToby Isaac 
4678a820b8eSToby Isaac  test:
4688a820b8eSToby Isaac    suffix: 1_continuous_trimmed
4698a820b8eSToby Isaac    args: -dim 1 -continuous 1 -trimmed 1
4708a820b8eSToby Isaac 
4718a820b8eSToby Isaac  test:
4728a820b8eSToby Isaac    suffix: 2_simplex_discontinuous_trimmed
4738a820b8eSToby Isaac    args: -dim 2 -tensor 0 -continuous 0 -trimmed 1
4748a820b8eSToby Isaac 
4758a820b8eSToby Isaac  test:
4768a820b8eSToby Isaac    suffix: 2_simplex_continuous_trimmed
4778a820b8eSToby Isaac    args: -dim 2 -tensor 0 -continuous 1 -trimmed 1
4788a820b8eSToby Isaac 
4798a820b8eSToby Isaac  test:
4808a820b8eSToby Isaac    suffix: 2_tensor_discontinuous_trimmed
4818a820b8eSToby Isaac    args: -dim 2 -tensor 1 -continuous 0 -trimmed 1
4828a820b8eSToby Isaac 
4838a820b8eSToby Isaac  test:
4848a820b8eSToby Isaac    suffix: 2_tensor_continuous_trimmed
4858a820b8eSToby Isaac    args: -dim 2 -tensor 1 -continuous 1 -trimmed 1
4868a820b8eSToby Isaac 
4878a820b8eSToby Isaac  test:
4888a820b8eSToby Isaac    suffix: 3_simplex_discontinuous_trimmed
4898a820b8eSToby Isaac    args: -dim 3 -tensor 0 -continuous 0 -trimmed 1
4908a820b8eSToby Isaac 
4918a820b8eSToby Isaac  test:
4928a820b8eSToby Isaac    suffix: 3_simplex_continuous_trimmed
4938a820b8eSToby Isaac    args: -dim 3 -tensor 0 -continuous 1 -trimmed 1
4948a820b8eSToby Isaac 
4958a820b8eSToby Isaac  test:
4968a820b8eSToby Isaac    suffix: 3_tensor_discontinuous_trimmed
4978a820b8eSToby Isaac    args: -dim 3 -tensor 1 -continuous 0 -trimmed 1
4988a820b8eSToby Isaac 
4998a820b8eSToby Isaac  test:
5008a820b8eSToby Isaac    suffix: 3_tensor_continuous_trimmed
5018a820b8eSToby Isaac    args: -dim 3 -tensor 1 -continuous 1 -trimmed 1
5028a820b8eSToby Isaac 
5038a820b8eSToby Isaac  test:
5048a820b8eSToby Isaac    suffix: 3_wedge_discontinuous_trimmed
5058a820b8eSToby Isaac    args: -dim 3 -tensor 2 -continuous 0 -trimmed 1
5068a820b8eSToby Isaac 
5078a820b8eSToby Isaac  test:
5088a820b8eSToby Isaac    suffix: 3_wedge_continuous_trimmed
5098a820b8eSToby Isaac    args: -dim 3 -tensor 2 -continuous 1 -trimmed 1
5108a820b8eSToby Isaac 
5118a820b8eSToby Isaac TEST*/
512