1*c4762a1bSJed Brown static char help[] = "Spectral element access patterns with Plex\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscdmplex.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown typedef struct { 6*c4762a1bSJed Brown PetscInt dim; /* Topological problem dimension */ 7*c4762a1bSJed Brown PetscInt Nf; /* Number of fields */ 8*c4762a1bSJed Brown PetscInt *Nc; /* Number of components per field */ 9*c4762a1bSJed Brown PetscInt *k; /* Spectral order per field */ 10*c4762a1bSJed Brown } AppCtx; 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 13*c4762a1bSJed Brown { 14*c4762a1bSJed Brown PetscInt len; 15*c4762a1bSJed Brown PetscBool flg; 16*c4762a1bSJed Brown PetscErrorCode ierr; 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown PetscFunctionBeginUser; 19*c4762a1bSJed Brown options->dim = 2; 20*c4762a1bSJed Brown options->Nf = 0; 21*c4762a1bSJed Brown options->Nc = NULL; 22*c4762a1bSJed Brown options->k = NULL; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "SEM Problem Options", "DMPLEX");CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "Problem dimension", "ex6.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-num_fields", "The number of fields", "ex6.c", options->Nf, &options->Nf, NULL,0);CHKERRQ(ierr); 27*c4762a1bSJed Brown if (options->Nf) { 28*c4762a1bSJed Brown len = options->Nf; 29*c4762a1bSJed Brown ierr = PetscMalloc1(len, &options->Nc);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = PetscOptionsIntArray("-num_components", "The number of components per field", "ex6.c", options->Nc, &len, &flg);CHKERRQ(ierr); 31*c4762a1bSJed Brown if (flg && (len != options->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %d should be %d", len, options->Nf); 32*c4762a1bSJed Brown len = options->Nf; 33*c4762a1bSJed Brown ierr = PetscMalloc1(len, &options->k);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = PetscOptionsIntArray("-order", "The spectral order per field", "ex6.c", options->k, &len, &flg);CHKERRQ(ierr); 35*c4762a1bSJed Brown if (flg && (len != options->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of order array is %d should be %d", len, options->Nf); 36*c4762a1bSJed Brown } 37*c4762a1bSJed Brown ierr = PetscOptionsEnd(); 38*c4762a1bSJed Brown PetscFunctionReturn(0); 39*c4762a1bSJed Brown } 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown static PetscErrorCode LoadData2D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt clSize, Vec u, AppCtx *user) 42*c4762a1bSJed Brown { 43*c4762a1bSJed Brown PetscInt i, j, f, c; 44*c4762a1bSJed Brown PetscErrorCode ierr; 45*c4762a1bSJed Brown PetscScalar *closure; 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown PetscFunctionBeginUser; 48*c4762a1bSJed Brown ierr = PetscMalloc1(clSize,&closure);CHKERRQ(ierr); 49*c4762a1bSJed Brown for (j = 0; j < Nj; ++j) { 50*c4762a1bSJed Brown for (i = 0; i < Ni; ++i) { 51*c4762a1bSJed Brown PetscInt ki, kj, o = 0; 52*c4762a1bSJed Brown ierr = PetscArrayzero(closure,clSize);CHKERRQ(ierr); 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown for (f = 0; f < user->Nf; ++f) { 55*c4762a1bSJed Brown PetscInt ioff = i*user->k[f], joff = j*user->k[f]; 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown for (kj = 0; kj <= user->k[f]; ++kj) { 58*c4762a1bSJed Brown for (ki = 0; ki <= user->k[f]; ++ki) { 59*c4762a1bSJed Brown for (c = 0; c < user->Nc[f]; ++c) { 60*c4762a1bSJed Brown closure[o++] = ((kj + joff)*(Ni*user->k[f]+1) + ki + ioff)*user->Nc[f]+c; 61*c4762a1bSJed Brown } 62*c4762a1bSJed Brown } 63*c4762a1bSJed Brown } 64*c4762a1bSJed Brown } 65*c4762a1bSJed Brown ierr = DMPlexVecSetClosure(dm, NULL, u, j*Ni+i, closure, INSERT_VALUES);CHKERRQ(ierr); 66*c4762a1bSJed Brown } 67*c4762a1bSJed Brown } 68*c4762a1bSJed Brown ierr = PetscFree(closure);CHKERRQ(ierr); 69*c4762a1bSJed Brown PetscFunctionReturn(0); 70*c4762a1bSJed Brown } 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown static PetscErrorCode LoadData3D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt Nk, PetscInt clSize, Vec u, AppCtx *user) 73*c4762a1bSJed Brown { 74*c4762a1bSJed Brown PetscInt i, j, k, f, c; 75*c4762a1bSJed Brown PetscErrorCode ierr; 76*c4762a1bSJed Brown PetscScalar *closure; 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown PetscFunctionBeginUser; 79*c4762a1bSJed Brown ierr = PetscMalloc1(clSize,&closure);CHKERRQ(ierr); 80*c4762a1bSJed Brown for (k = 0; k < Nk; ++k) { 81*c4762a1bSJed Brown for (j = 0; j < Nj; ++j) { 82*c4762a1bSJed Brown for (i = 0; i < Ni; ++i) { 83*c4762a1bSJed Brown PetscInt ki, kj, kk, o = 0; 84*c4762a1bSJed Brown ierr = PetscArrayzero(closure,clSize);CHKERRQ(ierr); 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown for (f = 0; f < user->Nf; ++f) { 87*c4762a1bSJed Brown PetscInt ioff = i*user->k[f], joff = j*user->k[f], koff = k*user->k[f]; 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown for (kk = 0; kk <= user->k[f]; ++kk) { 90*c4762a1bSJed Brown for (kj = 0; kj <= user->k[f]; ++kj) { 91*c4762a1bSJed Brown for (ki = 0; ki <= user->k[f]; ++ki) { 92*c4762a1bSJed Brown for (c = 0; c < user->Nc[f]; ++c) { 93*c4762a1bSJed Brown closure[o++] = (((kk + koff)*(Nj*user->k[f]+1) + kj + joff)*(Ni*user->k[f]+1) + ki + ioff)*user->Nc[f]+c; 94*c4762a1bSJed Brown } 95*c4762a1bSJed Brown } 96*c4762a1bSJed Brown } 97*c4762a1bSJed Brown } 98*c4762a1bSJed Brown } 99*c4762a1bSJed Brown ierr = DMPlexVecSetClosure(dm, NULL, u, (k*Nj+j)*Ni+i, closure, INSERT_VALUES);CHKERRQ(ierr); 100*c4762a1bSJed Brown } 101*c4762a1bSJed Brown } 102*c4762a1bSJed Brown } 103*c4762a1bSJed Brown ierr = PetscFree(closure);CHKERRQ(ierr); 104*c4762a1bSJed Brown PetscFunctionReturn(0); 105*c4762a1bSJed Brown } 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown static PetscErrorCode CheckPoint(DM dm, Vec u, PetscInt point, AppCtx *user) 108*c4762a1bSJed Brown { 109*c4762a1bSJed Brown PetscSection s; 110*c4762a1bSJed Brown PetscScalar *a; 111*c4762a1bSJed Brown const PetscScalar *array; 112*c4762a1bSJed Brown PetscInt dof, d; 113*c4762a1bSJed Brown PetscErrorCode ierr; 114*c4762a1bSJed Brown 115*c4762a1bSJed Brown PetscFunctionBeginUser; 116*c4762a1bSJed Brown ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 117*c4762a1bSJed Brown ierr = VecGetArrayRead(u, &array);CHKERRQ(ierr); 118*c4762a1bSJed Brown ierr = DMPlexPointLocalRead(dm, point, array, &a);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = PetscSectionGetDof(s, point, &dof);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Point %D: ", point);CHKERRQ(ierr); 121*c4762a1bSJed Brown for (d = 0; d < dof; ++d) { 122*c4762a1bSJed Brown if (d > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");CHKERRQ(ierr);} 123*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double) PetscRealPart(a[d]));CHKERRQ(ierr); 124*c4762a1bSJed Brown } 125*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = VecRestoreArrayRead(u, &array);CHKERRQ(ierr); 127*c4762a1bSJed Brown PetscFunctionReturn(0); 128*c4762a1bSJed Brown } 129*c4762a1bSJed Brown 130*c4762a1bSJed Brown static PetscErrorCode ReadData2D(DM dm, Vec u, AppCtx *user) 131*c4762a1bSJed Brown { 132*c4762a1bSJed Brown PetscInt cStart, cEnd, cell; 133*c4762a1bSJed Brown PetscErrorCode ierr; 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown PetscFunctionBeginUser; 136*c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 137*c4762a1bSJed Brown for (cell = cStart; cell < cEnd; ++cell) { 138*c4762a1bSJed Brown PetscScalar *closure = NULL; 139*c4762a1bSJed Brown PetscInt closureSize, ki, kj, f, c, foff = 0; 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown ierr = DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure);CHKERRQ(ierr); 142*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D\n", cell);CHKERRQ(ierr); 143*c4762a1bSJed Brown for (f = 0; f < user->Nf; ++f) { 144*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " Field %D\n", f);CHKERRQ(ierr); 145*c4762a1bSJed Brown for (kj = user->k[f]; kj >= 0; --kj) { 146*c4762a1bSJed Brown for (ki = 0; ki <= user->k[f]; ++ki) { 147*c4762a1bSJed Brown if (ki > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, " ");CHKERRQ(ierr);} 148*c4762a1bSJed Brown for (c = 0; c < user->Nc[f]; ++c) { 149*c4762a1bSJed Brown if (c > 0) ierr = PetscPrintf(PETSC_COMM_SELF, ",");CHKERRQ(ierr); 150*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double) PetscRealPart(closure[(kj*(user->k[f]+1) + ki)*user->Nc[f]+c + foff]));CHKERRQ(ierr); 151*c4762a1bSJed Brown } 152*c4762a1bSJed Brown } 153*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 154*c4762a1bSJed Brown } 155*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\n\n");CHKERRQ(ierr); 156*c4762a1bSJed Brown foff += PetscSqr(user->k[f]+1); 157*c4762a1bSJed Brown } 158*c4762a1bSJed Brown ierr = DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\n\n");CHKERRQ(ierr); 160*c4762a1bSJed Brown } 161*c4762a1bSJed Brown PetscFunctionReturn(0); 162*c4762a1bSJed Brown } 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown static PetscErrorCode ReadData3D(DM dm, Vec u, AppCtx *user) 165*c4762a1bSJed Brown { 166*c4762a1bSJed Brown PetscInt cStart, cEnd, cell; 167*c4762a1bSJed Brown PetscErrorCode ierr; 168*c4762a1bSJed Brown 169*c4762a1bSJed Brown PetscFunctionBeginUser; 170*c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 171*c4762a1bSJed Brown for (cell = cStart; cell < cEnd; ++cell) { 172*c4762a1bSJed Brown PetscScalar *closure = NULL; 173*c4762a1bSJed Brown PetscInt closureSize, ki, kj, kk, f, c, foff = 0; 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown ierr = DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure);CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "Cell %D\n", cell);CHKERRQ(ierr); 177*c4762a1bSJed Brown for (f = 0; f < user->Nf; ++f) { 178*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " Field %D\n", f);CHKERRQ(ierr); 179*c4762a1bSJed Brown for (kk = user->k[f]; kk >= 0; --kk) { 180*c4762a1bSJed Brown for (kj = user->k[f]; kj >= 0; --kj) { 181*c4762a1bSJed Brown for (ki = 0; ki <= user->k[f]; ++ki) { 182*c4762a1bSJed Brown if (ki > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, " ");CHKERRQ(ierr);} 183*c4762a1bSJed Brown for (c = 0; c < user->Nc[f]; ++c) { 184*c4762a1bSJed Brown if (c > 0) ierr = PetscPrintf(PETSC_COMM_SELF, ",");CHKERRQ(ierr); 185*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double) PetscRealPart(closure[((kk*(user->k[f]+1) + kj)*(user->k[f]+1) + ki)*user->Nc[f]+c + foff]));CHKERRQ(ierr); 186*c4762a1bSJed Brown } 187*c4762a1bSJed Brown } 188*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 189*c4762a1bSJed Brown } 190*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 191*c4762a1bSJed Brown } 192*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\n\n");CHKERRQ(ierr); 193*c4762a1bSJed Brown foff += PetscSqr(user->k[f]+1); 194*c4762a1bSJed Brown } 195*c4762a1bSJed Brown ierr = DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure);CHKERRQ(ierr); 196*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\n\n");CHKERRQ(ierr); 197*c4762a1bSJed Brown } 198*c4762a1bSJed Brown PetscFunctionReturn(0); 199*c4762a1bSJed Brown } 200*c4762a1bSJed Brown 201*c4762a1bSJed Brown static PetscErrorCode SetSymmetries(DM dm, PetscSection s, AppCtx *user) 202*c4762a1bSJed Brown { 203*c4762a1bSJed Brown PetscInt f, o, i, j, k, c, d; 204*c4762a1bSJed Brown DMLabel depthLabel; 205*c4762a1bSJed Brown PetscErrorCode ierr; 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown PetscFunctionBegin; 208*c4762a1bSJed Brown ierr = DMGetLabel(dm,"depth",&depthLabel);CHKERRQ(ierr); 209*c4762a1bSJed Brown for (f = 0; f < user->Nf; f++) { 210*c4762a1bSJed Brown PetscSectionSym sym; 211*c4762a1bSJed Brown 212*c4762a1bSJed Brown if (user->k[f] < 3) continue; /* No symmetries needed for order < 3, because no cell, facet, edge or vertex has more than one node */ 213*c4762a1bSJed Brown ierr = PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)s),depthLabel,&sym);CHKERRQ(ierr); 214*c4762a1bSJed Brown 215*c4762a1bSJed Brown for (d = 0; d <= user->dim; d++) { 216*c4762a1bSJed Brown if (d == 1) { 217*c4762a1bSJed Brown PetscInt numDof = user->k[f] - 1; 218*c4762a1bSJed Brown PetscInt numComp = user->Nc[f]; 219*c4762a1bSJed Brown PetscInt minOrnt = -2; 220*c4762a1bSJed Brown PetscInt maxOrnt = 2; 221*c4762a1bSJed Brown PetscInt **perms; 222*c4762a1bSJed Brown 223*c4762a1bSJed Brown ierr = PetscCalloc1(maxOrnt - minOrnt,&perms);CHKERRQ(ierr); 224*c4762a1bSJed Brown for (o = minOrnt; o < maxOrnt; o++) { 225*c4762a1bSJed Brown PetscInt *perm; 226*c4762a1bSJed Brown 227*c4762a1bSJed Brown if (o == -1 || !o) { /* identity */ 228*c4762a1bSJed Brown perms[o - minOrnt] = NULL; 229*c4762a1bSJed Brown } else { 230*c4762a1bSJed Brown ierr = PetscMalloc1(numDof * numComp, &perm);CHKERRQ(ierr); 231*c4762a1bSJed Brown for (i = numDof - 1, k = 0; i >= 0; i--) { 232*c4762a1bSJed Brown for (j = 0; j < numComp; j++, k++) perm[k] = i * numComp + j; 233*c4762a1bSJed Brown } 234*c4762a1bSJed Brown perms[o - minOrnt] = perm; 235*c4762a1bSJed Brown } 236*c4762a1bSJed Brown } 237*c4762a1bSJed Brown ierr = PetscSectionSymLabelSetStratum(sym,d,numDof*numComp,minOrnt,maxOrnt,PETSC_OWN_POINTER,(const PetscInt **) perms,NULL);CHKERRQ(ierr); 238*c4762a1bSJed Brown } else if (d == 2) { 239*c4762a1bSJed Brown PetscInt perEdge = user->k[f] - 1; 240*c4762a1bSJed Brown PetscInt numDof = perEdge * perEdge; 241*c4762a1bSJed Brown PetscInt numComp = user->Nc[f]; 242*c4762a1bSJed Brown PetscInt minOrnt = -4; 243*c4762a1bSJed Brown PetscInt maxOrnt = 4; 244*c4762a1bSJed Brown PetscInt **perms; 245*c4762a1bSJed Brown 246*c4762a1bSJed Brown ierr = PetscCalloc1(maxOrnt-minOrnt,&perms);CHKERRQ(ierr); 247*c4762a1bSJed Brown for (o = minOrnt; o < maxOrnt; o++) { 248*c4762a1bSJed Brown PetscInt *perm; 249*c4762a1bSJed Brown 250*c4762a1bSJed Brown if (!o) continue; /* identity */ 251*c4762a1bSJed Brown ierr = PetscMalloc1(numDof * numComp, &perm);CHKERRQ(ierr); 252*c4762a1bSJed Brown /* We want to perm[k] to list which *localArray* position the *sectionArray* position k should go to for the given orientation*/ 253*c4762a1bSJed Brown switch (o) { 254*c4762a1bSJed Brown case 0: 255*c4762a1bSJed Brown break; /* identity */ 256*c4762a1bSJed Brown case -4: /* flip along (-1,-1)--( 1, 1), which swaps edges 0 and 3 and edges 1 and 2. This swaps the i and j variables */ 257*c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 258*c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 259*c4762a1bSJed Brown for (c = 0; c < numComp; c++) { 260*c4762a1bSJed Brown perm[k * numComp + c] = (perEdge * j + i) * numComp + c; 261*c4762a1bSJed Brown } 262*c4762a1bSJed Brown } 263*c4762a1bSJed Brown } 264*c4762a1bSJed Brown break; 265*c4762a1bSJed Brown case -3: /* flip along (-1, 0)--( 1, 0), which swaps edges 0 and 2. This reverses the i variable */ 266*c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 267*c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 268*c4762a1bSJed Brown for (c = 0; c < numComp; c++) { 269*c4762a1bSJed Brown perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + j) * numComp + c; 270*c4762a1bSJed Brown } 271*c4762a1bSJed Brown } 272*c4762a1bSJed Brown } 273*c4762a1bSJed Brown break; 274*c4762a1bSJed Brown case -2: /* flip along ( 1,-1)--(-1, 1), which swaps edges 0 and 1 and edges 2 and 3. This swaps the i and j variables and reverse both */ 275*c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 276*c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 277*c4762a1bSJed Brown for (c = 0; c < numComp; c++) { 278*c4762a1bSJed Brown perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + (perEdge - 1 - i)) * numComp + c; 279*c4762a1bSJed Brown } 280*c4762a1bSJed Brown } 281*c4762a1bSJed Brown } 282*c4762a1bSJed Brown break; 283*c4762a1bSJed Brown case -1: /* flip along ( 0,-1)--( 0, 1), which swaps edges 3 and 1. This reverses the j variable */ 284*c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 285*c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 286*c4762a1bSJed Brown for (c = 0; c < numComp; c++) { 287*c4762a1bSJed Brown perm[k * numComp + c] = (perEdge * i + (perEdge - 1 - j)) * numComp + c; 288*c4762a1bSJed Brown } 289*c4762a1bSJed Brown } 290*c4762a1bSJed Brown } 291*c4762a1bSJed Brown break; 292*c4762a1bSJed Brown case 1: /* rotate section edge 1 to local edge 0. This swaps the i and j variables and then reverses the j variable */ 293*c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 294*c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 295*c4762a1bSJed Brown for (c = 0; c < numComp; c++) { 296*c4762a1bSJed Brown perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + i) * numComp + c; 297*c4762a1bSJed Brown } 298*c4762a1bSJed Brown } 299*c4762a1bSJed Brown } 300*c4762a1bSJed Brown break; 301*c4762a1bSJed Brown case 2: /* rotate section edge 2 to local edge 0. This reverse both i and j variables */ 302*c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 303*c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 304*c4762a1bSJed Brown for (c = 0; c < numComp; c++) { 305*c4762a1bSJed Brown perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + (perEdge - 1 - j)) * numComp + c; 306*c4762a1bSJed Brown } 307*c4762a1bSJed Brown } 308*c4762a1bSJed Brown } 309*c4762a1bSJed Brown break; 310*c4762a1bSJed Brown case 3: /* rotate section edge 3 to local edge 0. This swaps the i and j variables and then reverses the i variable */ 311*c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 312*c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 313*c4762a1bSJed Brown for (c = 0; c < numComp; c++) { 314*c4762a1bSJed Brown perm[k * numComp + c] = (perEdge * j + (perEdge - 1 - i)) * numComp + c; 315*c4762a1bSJed Brown } 316*c4762a1bSJed Brown } 317*c4762a1bSJed Brown } 318*c4762a1bSJed Brown break; 319*c4762a1bSJed Brown default: 320*c4762a1bSJed Brown break; 321*c4762a1bSJed Brown } 322*c4762a1bSJed Brown perms[o - minOrnt] = perm; 323*c4762a1bSJed Brown } 324*c4762a1bSJed Brown ierr = PetscSectionSymLabelSetStratum(sym,d,numDof*numComp,minOrnt,maxOrnt,PETSC_OWN_POINTER,(const PetscInt **) perms,NULL);CHKERRQ(ierr); 325*c4762a1bSJed Brown } 326*c4762a1bSJed Brown } 327*c4762a1bSJed Brown ierr = PetscSectionSetFieldSym(s,f,sym);CHKERRQ(ierr); 328*c4762a1bSJed Brown ierr = PetscSectionSymDestroy(&sym);CHKERRQ(ierr); 329*c4762a1bSJed Brown } 330*c4762a1bSJed Brown ierr = PetscSectionViewFromOptions(s,NULL,"-section_with_sym_view");CHKERRQ(ierr); 331*c4762a1bSJed Brown PetscFunctionReturn(0); 332*c4762a1bSJed Brown } 333*c4762a1bSJed Brown 334*c4762a1bSJed Brown int main(int argc, char **argv) 335*c4762a1bSJed Brown { 336*c4762a1bSJed Brown DM dm; 337*c4762a1bSJed Brown PetscSection s; 338*c4762a1bSJed Brown Vec u; 339*c4762a1bSJed Brown AppCtx user; 340*c4762a1bSJed Brown PetscInt cells[3] = {2, 2, 2}; 341*c4762a1bSJed Brown PetscInt size = 0, f; 342*c4762a1bSJed Brown PetscErrorCode ierr; 343*c4762a1bSJed Brown 344*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 345*c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 346*c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, user.dim, PETSC_FALSE, cells, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr); 347*c4762a1bSJed Brown ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 348*c4762a1bSJed Brown ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 349*c4762a1bSJed Brown /* Create a section for SEM order k */ 350*c4762a1bSJed Brown { 351*c4762a1bSJed Brown PetscInt *numDof, d; 352*c4762a1bSJed Brown 353*c4762a1bSJed Brown ierr = PetscMalloc1(user.Nf*(user.dim+1), &numDof);CHKERRQ(ierr); 354*c4762a1bSJed Brown for (f = 0; f < user.Nf; ++f) { 355*c4762a1bSJed Brown for (d = 0; d <= user.dim; ++d) numDof[f*(user.dim+1)+d] = PetscPowInt(user.k[f]-1, d)*user.Nc[f]; 356*c4762a1bSJed Brown size += PetscPowInt(user.k[f]+1, d)*user.Nc[f]; 357*c4762a1bSJed Brown } 358*c4762a1bSJed Brown ierr = DMSetNumFields(dm, user.Nf);CHKERRQ(ierr); 359*c4762a1bSJed Brown ierr = DMPlexCreateSection(dm, NULL, user.Nc, numDof, 0, NULL, NULL, NULL, NULL, &s);CHKERRQ(ierr); 360*c4762a1bSJed Brown ierr = SetSymmetries(dm, s, &user);CHKERRQ(ierr); 361*c4762a1bSJed Brown ierr = PetscFree(numDof);CHKERRQ(ierr); 362*c4762a1bSJed Brown } 363*c4762a1bSJed Brown ierr = DMSetLocalSection(dm, s);CHKERRQ(ierr); 364*c4762a1bSJed Brown /* Create spectral ordering and load in data */ 365*c4762a1bSJed Brown ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);CHKERRQ(ierr); 366*c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &u);CHKERRQ(ierr); 367*c4762a1bSJed Brown switch (user.dim) { 368*c4762a1bSJed Brown case 2: ierr = LoadData2D(dm, 2, 2, size, u, &user);CHKERRQ(ierr);break; 369*c4762a1bSJed Brown case 3: ierr = LoadData3D(dm, 2, 2, 2, size, u, &user);CHKERRQ(ierr);break; 370*c4762a1bSJed Brown } 371*c4762a1bSJed Brown /* Remove ordering and check some values */ 372*c4762a1bSJed Brown ierr = PetscSectionSetClosurePermutation(s, (PetscObject) dm, NULL);CHKERRQ(ierr); 373*c4762a1bSJed Brown switch (user.dim) { 374*c4762a1bSJed Brown case 2: 375*c4762a1bSJed Brown ierr = CheckPoint(dm, u, 0, &user);CHKERRQ(ierr); 376*c4762a1bSJed Brown ierr = CheckPoint(dm, u, 13, &user);CHKERRQ(ierr); 377*c4762a1bSJed Brown ierr = CheckPoint(dm, u, 15, &user);CHKERRQ(ierr); 378*c4762a1bSJed Brown ierr = CheckPoint(dm, u, 19, &user);CHKERRQ(ierr); 379*c4762a1bSJed Brown break; 380*c4762a1bSJed Brown case 3: 381*c4762a1bSJed Brown ierr = CheckPoint(dm, u, 0, &user);CHKERRQ(ierr); 382*c4762a1bSJed Brown ierr = CheckPoint(dm, u, 13, &user);CHKERRQ(ierr); 383*c4762a1bSJed Brown ierr = CheckPoint(dm, u, 15, &user);CHKERRQ(ierr); 384*c4762a1bSJed Brown ierr = CheckPoint(dm, u, 19, &user);CHKERRQ(ierr); 385*c4762a1bSJed Brown break; 386*c4762a1bSJed Brown } 387*c4762a1bSJed Brown /* Recreate spectral ordering and read out data */ 388*c4762a1bSJed Brown ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s);CHKERRQ(ierr); 389*c4762a1bSJed Brown switch (user.dim) { 390*c4762a1bSJed Brown case 2: ierr = ReadData2D(dm, u, &user);CHKERRQ(ierr);break; 391*c4762a1bSJed Brown case 3: ierr = ReadData3D(dm, u, &user);CHKERRQ(ierr);break; 392*c4762a1bSJed Brown } 393*c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &u);CHKERRQ(ierr); 394*c4762a1bSJed Brown ierr = PetscSectionDestroy(&s);CHKERRQ(ierr); 395*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 396*c4762a1bSJed Brown ierr = PetscFree(user.Nc);CHKERRQ(ierr); 397*c4762a1bSJed Brown ierr = PetscFree(user.k);CHKERRQ(ierr); 398*c4762a1bSJed Brown ierr = PetscFinalize(); 399*c4762a1bSJed Brown return ierr; 400*c4762a1bSJed Brown } 401*c4762a1bSJed Brown 402*c4762a1bSJed Brown /*TEST 403*c4762a1bSJed Brown 404*c4762a1bSJed Brown # Spectral ordering 2D 0-5 405*c4762a1bSJed Brown test: 406*c4762a1bSJed Brown suffix: 0 407*c4762a1bSJed Brown args: -dim 2 -num_fields 1 -num_components 1 -order 2 408*c4762a1bSJed Brown test: 409*c4762a1bSJed Brown suffix: 1 410*c4762a1bSJed Brown args: -dim 2 -num_fields 1 -num_components 1 -order 3 411*c4762a1bSJed Brown test: 412*c4762a1bSJed Brown suffix: 2 413*c4762a1bSJed Brown args: -dim 2 -num_fields 1 -num_components 1 -order 5 414*c4762a1bSJed Brown test: 415*c4762a1bSJed Brown suffix: 3 416*c4762a1bSJed Brown args: -dim 2 -num_fields 1 -num_components 2 -order 2 417*c4762a1bSJed Brown test: 418*c4762a1bSJed Brown suffix: 4 419*c4762a1bSJed Brown args: -dim 2 -num_fields 2 -num_components 1,1 -order 2,2 420*c4762a1bSJed Brown test: 421*c4762a1bSJed Brown suffix: 5 422*c4762a1bSJed Brown args: -dim 2 -num_fields 2 -num_components 1,2 -order 2,3 423*c4762a1bSJed Brown # Spectral ordering 3D 6-11 424*c4762a1bSJed Brown test: 425*c4762a1bSJed Brown suffix: 6 426*c4762a1bSJed Brown args: -dim 3 -num_fields 1 -num_components 1 -order 2 427*c4762a1bSJed Brown test: 428*c4762a1bSJed Brown suffix: 7 429*c4762a1bSJed Brown args: -dim 3 -num_fields 1 -num_components 1 -order 3 430*c4762a1bSJed Brown test: 431*c4762a1bSJed Brown suffix: 8 432*c4762a1bSJed Brown args: -dim 3 -num_fields 1 -num_components 1 -order 5 433*c4762a1bSJed Brown test: 434*c4762a1bSJed Brown suffix: 9 435*c4762a1bSJed Brown args: -dim 3 -num_fields 1 -num_components 2 -order 2 436*c4762a1bSJed Brown test: 437*c4762a1bSJed Brown suffix: 10 438*c4762a1bSJed Brown args: -dim 3 -num_fields 2 -num_components 1,1 -order 2,2 439*c4762a1bSJed Brown test: 440*c4762a1bSJed Brown suffix: 11 441*c4762a1bSJed Brown args: -dim 3 -num_fields 2 -num_components 1,2 -order 2,3 442*c4762a1bSJed Brown 443*c4762a1bSJed Brown TEST*/ 444