1*8a820b8eSToby Isaac 2*8a820b8eSToby Isaac #include <petscfe.h> 3*8a820b8eSToby Isaac #include <petscdmplex.h> 4*8a820b8eSToby Isaac #include <petsc/private/hashmap.h> 5*8a820b8eSToby Isaac #include <petsc/private/dmpleximpl.h> 6*8a820b8eSToby Isaac #include <petsc/private/petscfeimpl.h> 7*8a820b8eSToby Isaac 8*8a820b8eSToby Isaac const char help[] = "Test PETSCDUALSPACELAGRANGE\n"; 9*8a820b8eSToby Isaac 10*8a820b8eSToby Isaac typedef struct _PetscHashLagKey 11*8a820b8eSToby Isaac { 12*8a820b8eSToby Isaac PetscInt dim; 13*8a820b8eSToby Isaac PetscInt order; 14*8a820b8eSToby Isaac PetscInt formDegree; 15*8a820b8eSToby Isaac PetscBool trimmed; 16*8a820b8eSToby Isaac PetscBool tensor; 17*8a820b8eSToby Isaac PetscBool continuous; 18*8a820b8eSToby Isaac } PetscHashLagKey; 19*8a820b8eSToby Isaac 20*8a820b8eSToby Isaac #define PetscHashLagKeyHash(key) \ 21*8a820b8eSToby Isaac PetscHashCombine(PetscHashCombine(PetscHashCombine(PetscHashInt((key).dim), \ 22*8a820b8eSToby Isaac PetscHashInt((key).order)), \ 23*8a820b8eSToby Isaac PetscHashInt((key).formDegree)), \ 24*8a820b8eSToby Isaac PetscHashCombine(PetscHashCombine(PetscHashInt((key).trimmed), \ 25*8a820b8eSToby Isaac PetscHashInt((key).tensor)), \ 26*8a820b8eSToby Isaac PetscHashInt((key).continuous))) 27*8a820b8eSToby Isaac 28*8a820b8eSToby Isaac #define PetscHashLagKeyEqual(k1,k2) \ 29*8a820b8eSToby Isaac (((k1).dim == (k2).dim) ? \ 30*8a820b8eSToby Isaac ((k1).order == (k2).order) ? \ 31*8a820b8eSToby Isaac ((k1).formDegree == (k2).formDegree) ? \ 32*8a820b8eSToby Isaac ((k1).trimmed == (k2).trimmed) ? \ 33*8a820b8eSToby Isaac ((k1).tensor == (k2).tensor) ? \ 34*8a820b8eSToby Isaac ((k1).continuous == (k2).continuous) : 0 : 0 : 0 : 0 : 0) 35*8a820b8eSToby Isaac 36*8a820b8eSToby Isaac PETSC_HASH_MAP(HashLag, PetscHashLagKey, PetscInt, PetscHashLagKeyHash, PetscHashLagKeyEqual, 0) 37*8a820b8eSToby Isaac 38*8a820b8eSToby Isaac static PetscErrorCode ExpectedNumDofs_Total(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs); 39*8a820b8eSToby Isaac static PetscErrorCode ExpectedNumDofs_Interior(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs); 40*8a820b8eSToby Isaac 41*8a820b8eSToby Isaac static PetscErrorCode ExpectedNumDofs_Total(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscInt nCopies, PetscInt *nDofs) 42*8a820b8eSToby Isaac { 43*8a820b8eSToby Isaac PetscErrorCode ierr; 44*8a820b8eSToby Isaac 45*8a820b8eSToby Isaac PetscFunctionBegin; 46*8a820b8eSToby Isaac formDegree = PetscAbsInt(formDegree); 47*8a820b8eSToby Isaac /* see femtable.org for the source of most of these values */ 48*8a820b8eSToby Isaac *nDofs = -1; 49*8a820b8eSToby Isaac if (tensor == 0) { /* simplex spaces */ 50*8a820b8eSToby Isaac if (trimmed) { 51*8a820b8eSToby Isaac PetscInt rnchooserk; 52*8a820b8eSToby Isaac PetscInt rkm1choosek; 53*8a820b8eSToby Isaac 54*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(order + dim, order + formDegree, &rnchooserk);CHKERRQ(ierr); 55*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(order + formDegree - 1, formDegree, &rkm1choosek);CHKERRQ(ierr); 56*8a820b8eSToby Isaac *nDofs = rnchooserk * rkm1choosek * nCopies; 57*8a820b8eSToby Isaac } else { 58*8a820b8eSToby Isaac PetscInt rnchooserk; 59*8a820b8eSToby Isaac PetscInt rkchoosek; 60*8a820b8eSToby Isaac 61*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(order + dim, order + formDegree, &rnchooserk);CHKERRQ(ierr); 62*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(order + formDegree, formDegree, &rkchoosek);CHKERRQ(ierr); 63*8a820b8eSToby Isaac *nDofs = rnchooserk * rkchoosek * nCopies; 64*8a820b8eSToby Isaac } 65*8a820b8eSToby Isaac } else if (tensor == 1) { /* hypercubes */ 66*8a820b8eSToby Isaac if (trimmed) { 67*8a820b8eSToby Isaac PetscInt nchoosek; 68*8a820b8eSToby Isaac PetscInt rpowk, rp1pownmk; 69*8a820b8eSToby Isaac 70*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(dim, formDegree, &nchoosek);CHKERRQ(ierr); 71*8a820b8eSToby Isaac rpowk = PetscPowInt(order, formDegree);CHKERRQ(ierr); 72*8a820b8eSToby Isaac rp1pownmk = PetscPowInt(order + 1, dim - formDegree);CHKERRQ(ierr); 73*8a820b8eSToby Isaac *nDofs = nchoosek * rpowk * rp1pownmk * nCopies; 74*8a820b8eSToby Isaac } else { 75*8a820b8eSToby Isaac PetscInt nchoosek; 76*8a820b8eSToby Isaac PetscInt rp1pown; 77*8a820b8eSToby Isaac 78*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(dim, formDegree, &nchoosek);CHKERRQ(ierr); 79*8a820b8eSToby Isaac rp1pown = PetscPowInt(order + 1, dim);CHKERRQ(ierr); 80*8a820b8eSToby Isaac *nDofs = nchoosek * rp1pown * nCopies; 81*8a820b8eSToby Isaac } 82*8a820b8eSToby Isaac } else { /* prism */ 83*8a820b8eSToby Isaac PetscInt tracek = 0; 84*8a820b8eSToby Isaac PetscInt tracekm1 = 0; 85*8a820b8eSToby Isaac PetscInt fiber0 = 0; 86*8a820b8eSToby Isaac PetscInt fiber1 = 0; 87*8a820b8eSToby Isaac 88*8a820b8eSToby Isaac if (formDegree < dim) { 89*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Total(dim - 1, order, formDegree, trimmed, 0, 1, &tracek);CHKERRQ(ierr); 90*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Total(1, order, 0, trimmed, 0, 1, &fiber0);CHKERRQ(ierr); 91*8a820b8eSToby Isaac } 92*8a820b8eSToby Isaac if (formDegree > 0) { 93*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Total(dim - 1, order, formDegree - 1, trimmed, 0, 1, &tracekm1);CHKERRQ(ierr); 94*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Total(1, order, 1, trimmed, 0, 1, &fiber1);CHKERRQ(ierr); 95*8a820b8eSToby Isaac } 96*8a820b8eSToby Isaac *nDofs = (tracek * fiber0 + tracekm1 * fiber1) * nCopies; 97*8a820b8eSToby Isaac } 98*8a820b8eSToby Isaac PetscFunctionReturn(0); 99*8a820b8eSToby Isaac } 100*8a820b8eSToby Isaac 101*8a820b8eSToby Isaac static PetscErrorCode ExpectedNumDofs_Interior(PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, 102*8a820b8eSToby Isaac PetscInt tensor, PetscInt nCopies, PetscInt *nDofs) 103*8a820b8eSToby Isaac { 104*8a820b8eSToby Isaac PetscErrorCode ierr; 105*8a820b8eSToby Isaac 106*8a820b8eSToby Isaac PetscFunctionBegin; 107*8a820b8eSToby Isaac formDegree = PetscAbsInt(formDegree); 108*8a820b8eSToby Isaac /* see femtable.org for the source of most of these values */ 109*8a820b8eSToby Isaac *nDofs = -1; 110*8a820b8eSToby Isaac if (tensor == 0) { /* simplex spaces */ 111*8a820b8eSToby Isaac if (trimmed) { 112*8a820b8eSToby Isaac if (order + formDegree > dim) { 113*8a820b8eSToby Isaac PetscInt eorder = order + formDegree - dim - 1; 114*8a820b8eSToby Isaac PetscInt eformDegree = dim - formDegree; 115*8a820b8eSToby Isaac PetscInt rnchooserk; 116*8a820b8eSToby Isaac PetscInt rkchoosek; 117*8a820b8eSToby Isaac 118*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(eorder + dim, eorder + eformDegree, &rnchooserk);CHKERRQ(ierr); 119*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(eorder + eformDegree, eformDegree, &rkchoosek);CHKERRQ(ierr); 120*8a820b8eSToby Isaac *nDofs = rnchooserk * rkchoosek * nCopies; 121*8a820b8eSToby Isaac } else { 122*8a820b8eSToby Isaac *nDofs = 0; 123*8a820b8eSToby Isaac } 124*8a820b8eSToby Isaac 125*8a820b8eSToby Isaac } else { 126*8a820b8eSToby Isaac if (order + formDegree > dim) { 127*8a820b8eSToby Isaac PetscInt eorder = order + formDegree - dim; 128*8a820b8eSToby Isaac PetscInt eformDegree = dim - formDegree; 129*8a820b8eSToby Isaac PetscInt rnchooserk; 130*8a820b8eSToby Isaac PetscInt rkm1choosek; 131*8a820b8eSToby Isaac 132*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(eorder + dim, eorder + eformDegree, &rnchooserk);CHKERRQ(ierr); 133*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(eorder + eformDegree - 1, eformDegree, &rkm1choosek);CHKERRQ(ierr); 134*8a820b8eSToby Isaac *nDofs = rnchooserk * rkm1choosek * nCopies; 135*8a820b8eSToby Isaac } else { 136*8a820b8eSToby Isaac *nDofs = 0; 137*8a820b8eSToby Isaac } 138*8a820b8eSToby Isaac } 139*8a820b8eSToby Isaac } else if (tensor == 1) { /* hypercubes */ 140*8a820b8eSToby Isaac if (dim < 2) { 141*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(dim, order, formDegree, trimmed, 0, nCopies, nDofs);CHKERRQ(ierr); 142*8a820b8eSToby Isaac } else { 143*8a820b8eSToby Isaac PetscInt tracek = 0; 144*8a820b8eSToby Isaac PetscInt tracekm1 = 0; 145*8a820b8eSToby Isaac PetscInt fiber0 = 0; 146*8a820b8eSToby Isaac PetscInt fiber1 = 0; 147*8a820b8eSToby Isaac 148*8a820b8eSToby Isaac if (formDegree < dim) { 149*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(dim - 1, order, formDegree, trimmed, dim > 2 ? 1 : 0, 1, &tracek);CHKERRQ(ierr); 150*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(1, order, 0, trimmed, 0, 1, &fiber0);CHKERRQ(ierr); 151*8a820b8eSToby Isaac } 152*8a820b8eSToby Isaac if (formDegree > 0) { 153*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(dim - 1, order, formDegree - 1, trimmed, dim > 2 ? 1 : 0, 1, &tracekm1);CHKERRQ(ierr); 154*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(1, order, 1, trimmed, 0, 1, &fiber1);CHKERRQ(ierr); 155*8a820b8eSToby Isaac } 156*8a820b8eSToby Isaac *nDofs = (tracek * fiber0 + tracekm1 * fiber1) * nCopies; 157*8a820b8eSToby Isaac } 158*8a820b8eSToby Isaac } else { /* prism */ 159*8a820b8eSToby Isaac PetscInt tracek = 0; 160*8a820b8eSToby Isaac PetscInt tracekm1 = 0; 161*8a820b8eSToby Isaac PetscInt fiber0 = 0; 162*8a820b8eSToby Isaac PetscInt fiber1 = 0; 163*8a820b8eSToby Isaac 164*8a820b8eSToby Isaac if (formDegree < dim) { 165*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(dim - 1, order, formDegree, trimmed, 0, 1, &tracek);CHKERRQ(ierr); 166*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(1, order, 0, trimmed, 0, 1, &fiber0);CHKERRQ(ierr); 167*8a820b8eSToby Isaac } 168*8a820b8eSToby Isaac if (formDegree > 0) { 169*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(dim - 1, order, formDegree - 1, trimmed, 0, 1, &tracekm1);CHKERRQ(ierr); 170*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(1, order, 1, trimmed, 0, 1, &fiber1);CHKERRQ(ierr); 171*8a820b8eSToby Isaac } 172*8a820b8eSToby Isaac *nDofs = (tracek * fiber0 + tracekm1 * fiber1) * nCopies; 173*8a820b8eSToby Isaac } 174*8a820b8eSToby Isaac PetscFunctionReturn(0); 175*8a820b8eSToby Isaac } 176*8a820b8eSToby Isaac 177*8a820b8eSToby Isaac PetscErrorCode testLagrange(PetscHashLag lagTable, DM K, PetscInt dim, PetscInt order, PetscInt formDegree, PetscBool trimmed, PetscInt tensor, PetscBool continuous, PetscInt nCopies) 178*8a820b8eSToby Isaac { 179*8a820b8eSToby Isaac PetscDualSpace sp; 180*8a820b8eSToby Isaac MPI_Comm comm = PETSC_COMM_SELF; 181*8a820b8eSToby Isaac PetscInt Nk; 182*8a820b8eSToby Isaac PetscHashLagKey key; 183*8a820b8eSToby Isaac PetscHashIter iter; 184*8a820b8eSToby Isaac PetscBool missing; 185*8a820b8eSToby Isaac PetscInt spdim, spintdim, exspdim, exspintdim; 186*8a820b8eSToby Isaac PetscErrorCode ierr; 187*8a820b8eSToby Isaac 188*8a820b8eSToby Isaac PetscFunctionBegin; 189*8a820b8eSToby Isaac ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk);CHKERRQ(ierr); 190*8a820b8eSToby Isaac ierr = PetscDualSpaceCreate(comm, &sp);CHKERRQ(ierr); 191*8a820b8eSToby Isaac ierr = PetscDualSpaceSetType(sp, PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 192*8a820b8eSToby Isaac ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr); 193*8a820b8eSToby Isaac ierr = PetscDualSpaceSetOrder(sp, order);CHKERRQ(ierr); 194*8a820b8eSToby Isaac ierr = PetscDualSpaceSetFormDegree(sp, formDegree);CHKERRQ(ierr); 195*8a820b8eSToby Isaac ierr = PetscDualSpaceSetNumComponents(sp, nCopies * Nk);CHKERRQ(ierr); 196*8a820b8eSToby Isaac ierr = PetscDualSpaceLagrangeSetContinuity(sp, continuous);CHKERRQ(ierr); 197*8a820b8eSToby Isaac ierr = PetscDualSpaceLagrangeSetTensor(sp, (PetscBool) tensor);CHKERRQ(ierr); 198*8a820b8eSToby Isaac ierr = PetscDualSpaceLagrangeSetTrimmed(sp, trimmed);CHKERRQ(ierr); 199*8a820b8eSToby 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); 200*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Total(dim, order, formDegree, trimmed, tensor, nCopies, &exspdim);CHKERRQ(ierr); 201*8a820b8eSToby Isaac if (continuous && dim > 0 && order > 0) { 202*8a820b8eSToby Isaac ierr = ExpectedNumDofs_Interior(dim, order, formDegree, trimmed, tensor, nCopies, &exspintdim);CHKERRQ(ierr); 203*8a820b8eSToby Isaac } else { 204*8a820b8eSToby Isaac exspintdim = exspdim; 205*8a820b8eSToby Isaac } 206*8a820b8eSToby Isaac ierr = PetscDualSpaceSetUp(sp);CHKERRQ(ierr); 207*8a820b8eSToby Isaac ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr); 208*8a820b8eSToby Isaac ierr = PetscDualSpaceGetInteriorDimension(sp, &spintdim);CHKERRQ(ierr); 209*8a820b8eSToby Isaac if (spdim != exspdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected dual space dimension %D, got %D\n", exspdim, spdim);CHKERRQ(ierr); 210*8a820b8eSToby Isaac if (spintdim != exspintdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected dual space interior dimension %D, got %D\n", exspintdim, spintdim);CHKERRQ(ierr); 211*8a820b8eSToby Isaac key.dim = dim; 212*8a820b8eSToby Isaac key.formDegree = formDegree; 213*8a820b8eSToby Isaac ierr = PetscDualSpaceGetOrder(sp, &key.order);CHKERRQ(ierr); 214*8a820b8eSToby Isaac ierr = PetscDualSpaceLagrangeGetContinuity(sp, &key.continuous);CHKERRQ(ierr); 215*8a820b8eSToby Isaac if (tensor == 2) { 216*8a820b8eSToby Isaac key.tensor = 2; 217*8a820b8eSToby Isaac } else { 218*8a820b8eSToby Isaac ierr = PetscDualSpaceLagrangeGetTensor(sp, &key.tensor);CHKERRQ(ierr); 219*8a820b8eSToby Isaac } 220*8a820b8eSToby Isaac ierr = PetscDualSpaceLagrangeGetTrimmed(sp, &key.trimmed);CHKERRQ(ierr); 221*8a820b8eSToby 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); 222*8a820b8eSToby Isaac ierr = PetscHashLagPut(lagTable, key, &iter, &missing);CHKERRQ(ierr); 223*8a820b8eSToby Isaac if (missing) { 224*8a820b8eSToby Isaac DMPolytopeType type; 225*8a820b8eSToby Isaac 226*8a820b8eSToby Isaac ierr = DMPlexGetCellType(K, 0, &type);CHKERRQ(ierr); 227*8a820b8eSToby 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); 228*8a820b8eSToby Isaac ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 229*8a820b8eSToby Isaac { 230*8a820b8eSToby Isaac PetscQuadrature intNodes, allNodes; 231*8a820b8eSToby Isaac Mat intMat, allMat; 232*8a820b8eSToby Isaac MatInfo info; 233*8a820b8eSToby Isaac PetscInt i, j, nodeIdxDim, nodeVecDim, nNodes; 234*8a820b8eSToby Isaac const PetscInt *nodeIdx; 235*8a820b8eSToby Isaac const PetscReal *nodeVec; 236*8a820b8eSToby Isaac 237*8a820b8eSToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 238*8a820b8eSToby Isaac 239*8a820b8eSToby Isaac ierr = PetscLagNodeIndicesGetData_Internal(lag->allNodeIndices, &nodeIdxDim, &nodeVecDim, &nNodes, &nodeIdx, &nodeVec);CHKERRQ(ierr); 240*8a820b8eSToby Isaac if (nodeVecDim != Nk) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect nodeVecDim");CHKERRQ(ierr); 241*8a820b8eSToby Isaac if (nNodes != spdim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect nNodes");CHKERRQ(ierr); 242*8a820b8eSToby Isaac 243*8a820b8eSToby Isaac ierr = PetscDualSpaceGetAllData(sp, &allNodes, &allMat);CHKERRQ(ierr); 244*8a820b8eSToby Isaac 245*8a820b8eSToby Isaac ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "All nodes:\n");CHKERRQ(ierr); 246*8a820b8eSToby Isaac ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 247*8a820b8eSToby Isaac ierr = PetscQuadratureView(allNodes, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 248*8a820b8eSToby Isaac ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 249*8a820b8eSToby Isaac ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "All node indices:\n");CHKERRQ(ierr); 250*8a820b8eSToby Isaac for (i = 0; i < spdim; i++) { 251*8a820b8eSToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "(");CHKERRQ(ierr); 252*8a820b8eSToby Isaac for (j = 0; j < nodeIdxDim; j++) { 253*8a820b8eSToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " %D,", nodeIdx[i * nodeIdxDim + j]);CHKERRQ(ierr); 254*8a820b8eSToby Isaac } 255*8a820b8eSToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "): [");CHKERRQ(ierr); 256*8a820b8eSToby Isaac for (j = 0; j < nodeVecDim; j++) { 257*8a820b8eSToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " %g,", nodeVec[i * nodeVecDim + j]);CHKERRQ(ierr); 258*8a820b8eSToby Isaac } 259*8a820b8eSToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "]\n");CHKERRQ(ierr); 260*8a820b8eSToby Isaac } 261*8a820b8eSToby Isaac 262*8a820b8eSToby Isaac ierr = MatGetInfo(allMat, MAT_LOCAL, &info);CHKERRQ(ierr); 263*8a820b8eSToby Isaac ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "All matrix: %D nonzeros\n", (PetscInt) info.nz_used);CHKERRQ(ierr); 264*8a820b8eSToby Isaac 265*8a820b8eSToby Isaac ierr = PetscDualSpaceGetInteriorData(sp, &intNodes, &intMat);CHKERRQ(ierr); 266*8a820b8eSToby Isaac if (intMat && intMat != allMat) { 267*8a820b8eSToby Isaac PetscInt intNodeIdxDim, intNodeVecDim, intNnodes; 268*8a820b8eSToby Isaac const PetscInt *intNodeIdx; 269*8a820b8eSToby Isaac const PetscReal *intNodeVec; 270*8a820b8eSToby Isaac PetscBool same; 271*8a820b8eSToby Isaac 272*8a820b8eSToby Isaac ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior nodes:\n");CHKERRQ(ierr); 273*8a820b8eSToby Isaac ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 274*8a820b8eSToby Isaac ierr = PetscQuadratureView(intNodes, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 275*8a820b8eSToby Isaac ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 276*8a820b8eSToby Isaac 277*8a820b8eSToby Isaac ierr = MatGetInfo(intMat, MAT_LOCAL, &info);CHKERRQ(ierr); 278*8a820b8eSToby Isaac ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior matrix: %D nonzeros\n", (PetscInt) info.nz_used);CHKERRQ(ierr); 279*8a820b8eSToby Isaac ierr = PetscLagNodeIndicesGetData_Internal(lag->intNodeIndices, &intNodeIdxDim, &intNodeVecDim, &intNnodes, &intNodeIdx, &intNodeVec);CHKERRQ(ierr); 280*8a820b8eSToby Isaac if (intNodeIdxDim != nodeIdxDim || intNodeVecDim != nodeVecDim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node indices not the same shale as all node indices"); 281*8a820b8eSToby Isaac if (intNnodes != spintdim) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Incorrect interior nNodes");CHKERRQ(ierr); 282*8a820b8eSToby Isaac ierr = PetscArraycmp(intNodeIdx, nodeIdx, nodeIdxDim * intNnodes, &same);CHKERRQ(ierr); 283*8a820b8eSToby Isaac if (!same) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node indices not the same as start of all node indices"); 284*8a820b8eSToby Isaac ierr = PetscArraycmp(intNodeVec, nodeVec, nodeVecDim * intNnodes, &same);CHKERRQ(ierr); 285*8a820b8eSToby Isaac if (!same) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node vectors not the same as start of all node vectors"); 286*8a820b8eSToby Isaac } else if (intMat) { 287*8a820b8eSToby Isaac ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior data is the same as all data\n");CHKERRQ(ierr); 288*8a820b8eSToby Isaac if (intNodes != allNodes) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior nodes should be the same as all nodes"); 289*8a820b8eSToby Isaac if (lag->intNodeIndices != lag->allNodeIndices) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Interior node indices should be the same as all node indices"); 290*8a820b8eSToby Isaac } 291*8a820b8eSToby Isaac } 292*8a820b8eSToby Isaac if (dim <= 2 && spintdim) { 293*8a820b8eSToby Isaac PetscInt coneSize, o; 294*8a820b8eSToby Isaac 295*8a820b8eSToby Isaac ierr = DMPlexGetConeSize(K, 0, &coneSize);CHKERRQ(ierr); 296*8a820b8eSToby Isaac for (o = -coneSize; o < coneSize; o++) { 297*8a820b8eSToby Isaac Mat symMat; 298*8a820b8eSToby Isaac 299*8a820b8eSToby Isaac ierr = PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(sp, o, &symMat);CHKERRQ(ierr); 300*8a820b8eSToby Isaac ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF, "Interior node symmetry matrix for orientation %D:\n", o);CHKERRQ(ierr); 301*8a820b8eSToby Isaac ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 302*8a820b8eSToby Isaac ierr = MatView(symMat, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 303*8a820b8eSToby Isaac ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 304*8a820b8eSToby Isaac ierr = MatDestroy(&symMat);CHKERRQ(ierr); 305*8a820b8eSToby Isaac } 306*8a820b8eSToby Isaac } 307*8a820b8eSToby Isaac ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 308*8a820b8eSToby Isaac } 309*8a820b8eSToby Isaac ierr = PetscDualSpaceDestroy(&sp);CHKERRQ(ierr); 310*8a820b8eSToby Isaac PetscFunctionReturn(0); 311*8a820b8eSToby Isaac } 312*8a820b8eSToby Isaac 313*8a820b8eSToby Isaac static PetscErrorCode DMPlexCreateReferenceWedge(MPI_Comm comm, DM *refdm) 314*8a820b8eSToby Isaac { 315*8a820b8eSToby Isaac PetscInt dim = 3; 316*8a820b8eSToby Isaac DM rdm; 317*8a820b8eSToby Isaac PetscErrorCode ierr; 318*8a820b8eSToby Isaac 319*8a820b8eSToby Isaac PetscFunctionBeginUser; 320*8a820b8eSToby Isaac ierr = DMCreate(comm, &rdm);CHKERRQ(ierr); 321*8a820b8eSToby Isaac ierr = DMSetType(rdm, DMPLEX);CHKERRQ(ierr); 322*8a820b8eSToby Isaac ierr = DMSetDimension(rdm, dim);CHKERRQ(ierr); 323*8a820b8eSToby Isaac { 324*8a820b8eSToby Isaac PetscInt numPoints[4] = {6, 9, 5, 1}; 325*8a820b8eSToby Isaac PetscInt coneSize[21] = {5, 326*8a820b8eSToby Isaac 3, 3, 327*8a820b8eSToby Isaac 4, 4, 4, 328*8a820b8eSToby Isaac 2, 2, 2, 2, 2, 2, 2, 2, 2, 329*8a820b8eSToby Isaac 0, 0, 0, 0, 0, 0}; 330*8a820b8eSToby Isaac PetscInt cones[41] = {1, 2, 3, 4, 5, 331*8a820b8eSToby Isaac 6, 7, 8, 332*8a820b8eSToby Isaac 9, 10, 11, 333*8a820b8eSToby Isaac 8, 12, 9, 13, 334*8a820b8eSToby Isaac 7, 14, 10, 12, 335*8a820b8eSToby Isaac 6, 13, 11, 14, 336*8a820b8eSToby Isaac 15, 16, 16, 17, 17, 15, 337*8a820b8eSToby Isaac 18, 19, 19, 20, 20, 18, 338*8a820b8eSToby Isaac 17, 19, 18, 15, 16, 20}; 339*8a820b8eSToby Isaac PetscInt coneOrientations[41] = {0, 0, 0, 0, 0, 340*8a820b8eSToby Isaac 0, 0, 0, 341*8a820b8eSToby Isaac 0, 0, 0, 342*8a820b8eSToby Isaac -2, 0, -2, 0, 343*8a820b8eSToby Isaac -2, 0, -2, -2, 344*8a820b8eSToby Isaac -2, -2, -2, -2, 345*8a820b8eSToby Isaac 0, 0, 0, 0, 0, 0, 346*8a820b8eSToby Isaac 0, 0, 0, 0, 0, 0, 347*8a820b8eSToby Isaac 0, 0, 0, 0, 0, 0}; 348*8a820b8eSToby Isaac PetscScalar vertexCoords[18] = {-1.0, -1.0, -1.0, 349*8a820b8eSToby Isaac -1.0, 1.0, -1.0, 350*8a820b8eSToby Isaac 1.0, -1.0, -1.0, 351*8a820b8eSToby Isaac -1.0, -1.0, 1.0, 352*8a820b8eSToby Isaac 1.0, -1.0, 1.0, 353*8a820b8eSToby Isaac -1.0, 1.0, 1.0}; 354*8a820b8eSToby Isaac 355*8a820b8eSToby Isaac ierr = DMPlexCreateFromDAG(rdm, 3, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr); 356*8a820b8eSToby Isaac } 357*8a820b8eSToby Isaac *refdm = rdm; 358*8a820b8eSToby Isaac PetscFunctionReturn(0); 359*8a820b8eSToby Isaac } 360*8a820b8eSToby Isaac 361*8a820b8eSToby Isaac int main (int argc, char **argv) 362*8a820b8eSToby Isaac { 363*8a820b8eSToby Isaac PetscInt dim; 364*8a820b8eSToby Isaac PetscHashLag lagTable; 365*8a820b8eSToby Isaac PetscInt tensorCell; 366*8a820b8eSToby Isaac PetscInt order, ordermin, ordermax; 367*8a820b8eSToby Isaac PetscBool continuous; 368*8a820b8eSToby Isaac PetscBool trimmed; 369*8a820b8eSToby Isaac DM dm; 370*8a820b8eSToby Isaac PetscErrorCode ierr; 371*8a820b8eSToby Isaac 372*8a820b8eSToby Isaac ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 373*8a820b8eSToby Isaac dim = 3; 374*8a820b8eSToby Isaac tensorCell = 0; 375*8a820b8eSToby Isaac continuous = PETSC_FALSE; 376*8a820b8eSToby Isaac trimmed = PETSC_FALSE; 377*8a820b8eSToby Isaac ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for PETSCDUALSPACELAGRANGE test","none");CHKERRQ(ierr); 378*8a820b8eSToby Isaac ierr = PetscOptionsRangeInt("-dim", "The spatial dimension","ex1.c",dim,&dim,NULL,0,3);CHKERRQ(ierr); 379*8a820b8eSToby Isaac ierr = PetscOptionsRangeInt("-tensor", "(0) simplex (1) hypercube (2) wedge","ex1.c",tensorCell,&tensorCell,NULL,0,3);CHKERRQ(ierr); 380*8a820b8eSToby Isaac ierr = PetscOptionsBool("-continuous", "Whether the dual space has continuity","ex1.c",continuous,&continuous,NULL);CHKERRQ(ierr); 381*8a820b8eSToby Isaac ierr = PetscOptionsBool("-trimmed", "Whether the dual space matches a trimmed polynomial space","ex1.c",trimmed,&trimmed,NULL);CHKERRQ(ierr); 382*8a820b8eSToby Isaac ierr = PetscOptionsEnd(); 383*8a820b8eSToby Isaac ierr = PetscHashLagCreate(&lagTable);CHKERRQ(ierr); 384*8a820b8eSToby Isaac 385*8a820b8eSToby Isaac if (tensorCell < 2) { 386*8a820b8eSToby Isaac ierr = DMPlexCreateReferenceCell(PETSC_COMM_SELF, dim, (PetscBool) !tensorCell, &dm);CHKERRQ(ierr); 387*8a820b8eSToby Isaac } else { 388*8a820b8eSToby Isaac ierr = DMPlexCreateReferenceWedge(PETSC_COMM_SELF, &dm);CHKERRQ(ierr); 389*8a820b8eSToby Isaac } 390*8a820b8eSToby Isaac ordermin = trimmed ? 1 : 0; 391*8a820b8eSToby Isaac ordermax = tensorCell == 2 ? 4 : tensorCell == 1 ? 3 : dim + 2; 392*8a820b8eSToby Isaac for (order = ordermin; order <= ordermax; order++) { 393*8a820b8eSToby Isaac PetscInt formDegree; 394*8a820b8eSToby Isaac 395*8a820b8eSToby Isaac for (formDegree = PetscMin(0,-dim+1); formDegree <= dim; formDegree++) { 396*8a820b8eSToby Isaac PetscInt nCopies; 397*8a820b8eSToby Isaac 398*8a820b8eSToby Isaac for (nCopies = 1; nCopies <= 3; nCopies++) { 399*8a820b8eSToby Isaac ierr = testLagrange(lagTable, dm, dim, order, formDegree, trimmed, (PetscBool) tensorCell, continuous, nCopies);CHKERRQ(ierr); 400*8a820b8eSToby Isaac } 401*8a820b8eSToby Isaac } 402*8a820b8eSToby Isaac } 403*8a820b8eSToby Isaac ierr = DMDestroy(&dm);CHKERRQ(ierr); 404*8a820b8eSToby Isaac ierr = PetscHashLagDestroy(&lagTable);CHKERRQ(ierr); 405*8a820b8eSToby Isaac ierr = PetscFinalize(); 406*8a820b8eSToby Isaac return ierr; 407*8a820b8eSToby Isaac } 408*8a820b8eSToby Isaac 409*8a820b8eSToby Isaac /*TEST 410*8a820b8eSToby Isaac 411*8a820b8eSToby Isaac test: 412*8a820b8eSToby Isaac suffix: 0 413*8a820b8eSToby Isaac args: -dim 0 414*8a820b8eSToby Isaac 415*8a820b8eSToby Isaac test: 416*8a820b8eSToby Isaac suffix: 1_discontinuous_full 417*8a820b8eSToby Isaac args: -dim 1 -continuous 0 -trimmed 0 418*8a820b8eSToby Isaac 419*8a820b8eSToby Isaac test: 420*8a820b8eSToby Isaac suffix: 1_continuous_full 421*8a820b8eSToby Isaac args: -dim 1 -continuous 1 -trimmed 0 422*8a820b8eSToby Isaac 423*8a820b8eSToby Isaac test: 424*8a820b8eSToby Isaac suffix: 2_simplex_discontinuous_full 425*8a820b8eSToby Isaac args: -dim 2 -tensor 0 -continuous 0 -trimmed 0 426*8a820b8eSToby Isaac 427*8a820b8eSToby Isaac test: 428*8a820b8eSToby Isaac suffix: 2_simplex_continuous_full 429*8a820b8eSToby Isaac args: -dim 2 -tensor 0 -continuous 1 -trimmed 0 430*8a820b8eSToby Isaac 431*8a820b8eSToby Isaac test: 432*8a820b8eSToby Isaac suffix: 2_tensor_discontinuous_full 433*8a820b8eSToby Isaac args: -dim 2 -tensor 1 -continuous 0 -trimmed 0 434*8a820b8eSToby Isaac 435*8a820b8eSToby Isaac test: 436*8a820b8eSToby Isaac suffix: 2_tensor_continuous_full 437*8a820b8eSToby Isaac args: -dim 2 -tensor 1 -continuous 1 -trimmed 0 438*8a820b8eSToby Isaac 439*8a820b8eSToby Isaac test: 440*8a820b8eSToby Isaac suffix: 3_simplex_discontinuous_full 441*8a820b8eSToby Isaac args: -dim 3 -tensor 0 -continuous 0 -trimmed 0 442*8a820b8eSToby Isaac 443*8a820b8eSToby Isaac test: 444*8a820b8eSToby Isaac suffix: 3_simplex_continuous_full 445*8a820b8eSToby Isaac args: -dim 3 -tensor 0 -continuous 1 -trimmed 0 446*8a820b8eSToby Isaac 447*8a820b8eSToby Isaac test: 448*8a820b8eSToby Isaac suffix: 3_tensor_discontinuous_full 449*8a820b8eSToby Isaac args: -dim 3 -tensor 1 -continuous 0 -trimmed 0 450*8a820b8eSToby Isaac 451*8a820b8eSToby Isaac test: 452*8a820b8eSToby Isaac suffix: 3_tensor_continuous_full 453*8a820b8eSToby Isaac args: -dim 3 -tensor 1 -continuous 1 -trimmed 0 454*8a820b8eSToby Isaac 455*8a820b8eSToby Isaac test: 456*8a820b8eSToby Isaac suffix: 3_wedge_discontinuous_full 457*8a820b8eSToby Isaac args: -dim 3 -tensor 2 -continuous 0 -trimmed 0 458*8a820b8eSToby Isaac 459*8a820b8eSToby Isaac test: 460*8a820b8eSToby Isaac suffix: 3_wedge_continuous_full 461*8a820b8eSToby Isaac args: -dim 3 -tensor 2 -continuous 1 -trimmed 0 462*8a820b8eSToby Isaac 463*8a820b8eSToby Isaac test: 464*8a820b8eSToby Isaac suffix: 1_discontinuous_trimmed 465*8a820b8eSToby Isaac args: -dim 1 -continuous 0 -trimmed 1 466*8a820b8eSToby Isaac 467*8a820b8eSToby Isaac test: 468*8a820b8eSToby Isaac suffix: 1_continuous_trimmed 469*8a820b8eSToby Isaac args: -dim 1 -continuous 1 -trimmed 1 470*8a820b8eSToby Isaac 471*8a820b8eSToby Isaac test: 472*8a820b8eSToby Isaac suffix: 2_simplex_discontinuous_trimmed 473*8a820b8eSToby Isaac args: -dim 2 -tensor 0 -continuous 0 -trimmed 1 474*8a820b8eSToby Isaac 475*8a820b8eSToby Isaac test: 476*8a820b8eSToby Isaac suffix: 2_simplex_continuous_trimmed 477*8a820b8eSToby Isaac args: -dim 2 -tensor 0 -continuous 1 -trimmed 1 478*8a820b8eSToby Isaac 479*8a820b8eSToby Isaac test: 480*8a820b8eSToby Isaac suffix: 2_tensor_discontinuous_trimmed 481*8a820b8eSToby Isaac args: -dim 2 -tensor 1 -continuous 0 -trimmed 1 482*8a820b8eSToby Isaac 483*8a820b8eSToby Isaac test: 484*8a820b8eSToby Isaac suffix: 2_tensor_continuous_trimmed 485*8a820b8eSToby Isaac args: -dim 2 -tensor 1 -continuous 1 -trimmed 1 486*8a820b8eSToby Isaac 487*8a820b8eSToby Isaac test: 488*8a820b8eSToby Isaac suffix: 3_simplex_discontinuous_trimmed 489*8a820b8eSToby Isaac args: -dim 3 -tensor 0 -continuous 0 -trimmed 1 490*8a820b8eSToby Isaac 491*8a820b8eSToby Isaac test: 492*8a820b8eSToby Isaac suffix: 3_simplex_continuous_trimmed 493*8a820b8eSToby Isaac args: -dim 3 -tensor 0 -continuous 1 -trimmed 1 494*8a820b8eSToby Isaac 495*8a820b8eSToby Isaac test: 496*8a820b8eSToby Isaac suffix: 3_tensor_discontinuous_trimmed 497*8a820b8eSToby Isaac args: -dim 3 -tensor 1 -continuous 0 -trimmed 1 498*8a820b8eSToby Isaac 499*8a820b8eSToby Isaac test: 500*8a820b8eSToby Isaac suffix: 3_tensor_continuous_trimmed 501*8a820b8eSToby Isaac args: -dim 3 -tensor 1 -continuous 1 -trimmed 1 502*8a820b8eSToby Isaac 503*8a820b8eSToby Isaac test: 504*8a820b8eSToby Isaac suffix: 3_wedge_discontinuous_trimmed 505*8a820b8eSToby Isaac args: -dim 3 -tensor 2 -continuous 0 -trimmed 1 506*8a820b8eSToby Isaac 507*8a820b8eSToby Isaac test: 508*8a820b8eSToby Isaac suffix: 3_wedge_continuous_trimmed 509*8a820b8eSToby Isaac args: -dim 3 -tensor 2 -continuous 1 -trimmed 1 510*8a820b8eSToby Isaac 511*8a820b8eSToby Isaac TEST*/ 512