xref: /petsc/src/dm/impls/plex/tutorials/ex6.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Spectral element access patterns with Plex\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown typedef struct {
6c4762a1bSJed Brown   PetscInt  Nf;  /* Number of fields */
7c4762a1bSJed Brown   PetscInt *Nc;  /* Number of components per field */
8c4762a1bSJed Brown   PetscInt *k;   /* Spectral order per field */
9c4762a1bSJed Brown } AppCtx;
10c4762a1bSJed Brown 
11c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   PetscInt       len;
14c4762a1bSJed Brown   PetscBool      flg;
15c4762a1bSJed Brown   PetscErrorCode ierr;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   PetscFunctionBeginUser;
18c4762a1bSJed Brown   options->Nf = 0;
19c4762a1bSJed Brown   options->Nc = NULL;
20c4762a1bSJed Brown   options->k  = NULL;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "SEM Problem Options", "DMPLEX");CHKERRQ(ierr);
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBoundedInt("-num_fields", "The number of fields", "ex6.c", options->Nf, &options->Nf, NULL,0));
24c4762a1bSJed Brown   if (options->Nf) {
25c4762a1bSJed Brown     len  = options->Nf;
265f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(len, &options->Nc));
275f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsIntArray("-num_components", "The number of components per field", "ex6.c", options->Nc, &len, &flg));
282c71b3e2SJacob Faibussowitsch     PetscCheckFalse(flg && (len != options->Nf),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %d should be %d", len, options->Nf);
29c4762a1bSJed Brown     len  = options->Nf;
305f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(len, &options->k));
315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsIntArray("-order", "The spectral order per field", "ex6.c", options->k, &len, &flg));
322c71b3e2SJacob Faibussowitsch     PetscCheckFalse(flg && (len != options->Nf),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of order array is %d should be %d", len, options->Nf);
33c4762a1bSJed Brown   }
341e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
35c4762a1bSJed Brown   PetscFunctionReturn(0);
36c4762a1bSJed Brown }
37c4762a1bSJed Brown 
38c4762a1bSJed Brown static PetscErrorCode LoadData2D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt clSize, Vec u, AppCtx *user)
39c4762a1bSJed Brown {
40c4762a1bSJed Brown   PetscInt       i, j, f, c;
41c4762a1bSJed Brown   PetscScalar *closure;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   PetscFunctionBeginUser;
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(clSize,&closure));
45c4762a1bSJed Brown   for (j = 0; j < Nj; ++j) {
46c4762a1bSJed Brown     for (i = 0; i < Ni; ++i) {
47c4762a1bSJed Brown       PetscInt    ki, kj, o = 0;
485f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscArrayzero(closure,clSize));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown       for (f = 0; f < user->Nf; ++f) {
51c4762a1bSJed Brown         PetscInt ioff = i*user->k[f], joff = j*user->k[f];
52c4762a1bSJed Brown 
53c4762a1bSJed Brown         for (kj = 0; kj <= user->k[f]; ++kj) {
54c4762a1bSJed Brown           for (ki = 0; ki <= user->k[f]; ++ki) {
55c4762a1bSJed Brown             for (c = 0; c < user->Nc[f]; ++c) {
56c4762a1bSJed Brown               closure[o++] = ((kj + joff)*(Ni*user->k[f]+1) + ki + ioff)*user->Nc[f]+c;
57c4762a1bSJed Brown             }
58c4762a1bSJed Brown           }
59c4762a1bSJed Brown         }
60c4762a1bSJed Brown       }
615f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexVecSetClosure(dm, NULL, u, j*Ni+i, closure, INSERT_VALUES));
62c4762a1bSJed Brown     }
63c4762a1bSJed Brown   }
645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(closure));
65c4762a1bSJed Brown   PetscFunctionReturn(0);
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown static PetscErrorCode LoadData3D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt Nk, PetscInt clSize, Vec u, AppCtx *user)
69c4762a1bSJed Brown {
70c4762a1bSJed Brown   PetscInt       i, j, k, f, c;
71c4762a1bSJed Brown   PetscScalar *closure;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   PetscFunctionBeginUser;
745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(clSize,&closure));
75c4762a1bSJed Brown   for (k = 0; k < Nk; ++k) {
76c4762a1bSJed Brown     for (j = 0; j < Nj; ++j) {
77c4762a1bSJed Brown       for (i = 0; i < Ni; ++i) {
78c4762a1bSJed Brown         PetscInt    ki, kj, kk, o = 0;
795f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscArrayzero(closure,clSize));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown         for (f = 0; f < user->Nf; ++f) {
82c4762a1bSJed Brown           PetscInt ioff = i*user->k[f], joff = j*user->k[f], koff = k*user->k[f];
83c4762a1bSJed Brown 
84c4762a1bSJed Brown           for (kk = 0; kk <= user->k[f]; ++kk) {
85c4762a1bSJed Brown             for (kj = 0; kj <= user->k[f]; ++kj) {
86c4762a1bSJed Brown               for (ki = 0; ki <= user->k[f]; ++ki) {
87c4762a1bSJed Brown                 for (c = 0; c < user->Nc[f]; ++c) {
88c4762a1bSJed Brown                   closure[o++] = (((kk + koff)*(Nj*user->k[f]+1) + kj + joff)*(Ni*user->k[f]+1) + ki + ioff)*user->Nc[f]+c;
89c4762a1bSJed Brown                 }
90c4762a1bSJed Brown               }
91c4762a1bSJed Brown             }
92c4762a1bSJed Brown           }
93c4762a1bSJed Brown         }
945f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexVecSetClosure(dm, NULL, u, (k*Nj+j)*Ni+i, closure, INSERT_VALUES));
95c4762a1bSJed Brown       }
96c4762a1bSJed Brown     }
97c4762a1bSJed Brown   }
985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(closure));
99c4762a1bSJed Brown   PetscFunctionReturn(0);
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown static PetscErrorCode CheckPoint(DM dm, Vec u, PetscInt point, AppCtx *user)
103c4762a1bSJed Brown {
104c4762a1bSJed Brown   PetscSection       s;
105c4762a1bSJed Brown   PetscScalar        *a;
106c4762a1bSJed Brown   const PetscScalar  *array;
107c4762a1bSJed Brown   PetscInt           dof, d;
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   PetscFunctionBeginUser;
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dm, &s));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(u, &array));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexPointLocalRead(dm, point, array, &a));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetDof(s, point, &dof));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Point %D: ", point));
115c4762a1bSJed Brown   for (d = 0; d < dof; ++d) {
1165f80ce2aSJacob Faibussowitsch     if (d > 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, ", "));
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double) PetscRealPart(a[d])));
118c4762a1bSJed Brown   }
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n"));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(u, &array));
121c4762a1bSJed Brown   PetscFunctionReturn(0);
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown static PetscErrorCode ReadData2D(DM dm, Vec u, AppCtx *user)
125c4762a1bSJed Brown {
126c4762a1bSJed Brown   PetscInt       cStart, cEnd, cell;
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   PetscFunctionBeginUser;
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
130c4762a1bSJed Brown   for (cell = cStart; cell < cEnd; ++cell) {
131c4762a1bSJed Brown     PetscScalar *closure = NULL;
132c4762a1bSJed Brown     PetscInt     closureSize, ki, kj, f, c, foff = 0;
133c4762a1bSJed Brown 
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure));
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Cell %D\n", cell));
136c4762a1bSJed Brown     for (f = 0; f < user->Nf; ++f) {
1375f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "  Field %D\n", f));
138c4762a1bSJed Brown       for (kj = user->k[f]; kj >= 0; --kj) {
139c4762a1bSJed Brown         for (ki = 0; ki <= user->k[f]; ++ki) {
1405f80ce2aSJacob Faibussowitsch           if (ki > 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "  "));
141c4762a1bSJed Brown           for (c = 0; c < user->Nc[f]; ++c) {
1425f80ce2aSJacob Faibussowitsch             if (c > 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, ","));
1435f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double) PetscRealPart(closure[(kj*(user->k[f]+1) + ki)*user->Nc[f]+c + foff])));
144c4762a1bSJed Brown           }
145c4762a1bSJed Brown         }
1465f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n"));
147c4762a1bSJed Brown       }
1485f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
149c4762a1bSJed Brown       foff += PetscSqr(user->k[f]+1);
150c4762a1bSJed Brown     }
1515f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure));
1525f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
153c4762a1bSJed Brown   }
154c4762a1bSJed Brown   PetscFunctionReturn(0);
155c4762a1bSJed Brown }
156c4762a1bSJed Brown 
157c4762a1bSJed Brown static PetscErrorCode ReadData3D(DM dm, Vec u, AppCtx *user)
158c4762a1bSJed Brown {
159c4762a1bSJed Brown   PetscInt       cStart, cEnd, cell;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   PetscFunctionBeginUser;
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
163c4762a1bSJed Brown   for (cell = cStart; cell < cEnd; ++cell) {
164c4762a1bSJed Brown     PetscScalar *closure = NULL;
165c4762a1bSJed Brown     PetscInt     closureSize, ki, kj, kk, f, c, foff = 0;
166c4762a1bSJed Brown 
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure));
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Cell %D\n", cell));
169c4762a1bSJed Brown     for (f = 0; f < user->Nf; ++f) {
1705f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "  Field %D\n", f));
171c4762a1bSJed Brown       for (kk = user->k[f]; kk >= 0; --kk) {
172c4762a1bSJed Brown         for (kj = user->k[f]; kj >= 0; --kj) {
173c4762a1bSJed Brown           for (ki = 0; ki <= user->k[f]; ++ki) {
1745f80ce2aSJacob Faibussowitsch             if (ki > 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "  "));
175c4762a1bSJed Brown             for (c = 0; c < user->Nc[f]; ++c) {
1765f80ce2aSJacob Faibussowitsch               if (c > 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, ","));
1775f80ce2aSJacob Faibussowitsch               CHKERRQ(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])));
178c4762a1bSJed Brown             }
179c4762a1bSJed Brown           }
1805f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n"));
181c4762a1bSJed Brown         }
1825f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n"));
183c4762a1bSJed Brown       }
1845f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
185c4762a1bSJed Brown       foff += PetscSqr(user->k[f]+1);
186c4762a1bSJed Brown     }
1875f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure));
1885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
189c4762a1bSJed Brown   }
190c4762a1bSJed Brown   PetscFunctionReturn(0);
191c4762a1bSJed Brown }
192c4762a1bSJed Brown 
193c4762a1bSJed Brown static PetscErrorCode SetSymmetries(DM dm, PetscSection s, AppCtx *user)
194c4762a1bSJed Brown {
19530602db0SMatthew G. Knepley   PetscInt       dim, f, o, i, j, k, c, d;
196c4762a1bSJed Brown   DMLabel        depthLabel;
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   PetscFunctionBegin;
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm,"depth",&depthLabel));
201c4762a1bSJed Brown   for (f = 0; f < user->Nf; f++) {
202c4762a1bSJed Brown     PetscSectionSym sym;
203c4762a1bSJed Brown 
204c4762a1bSJed 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 */
2055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)s),depthLabel,&sym));
206c4762a1bSJed Brown 
20730602db0SMatthew G. Knepley     for (d = 0; d <= dim; d++) {
208c4762a1bSJed Brown       if (d == 1) {
209c4762a1bSJed Brown         PetscInt        numDof  = user->k[f] - 1;
210c4762a1bSJed Brown         PetscInt        numComp = user->Nc[f];
211b5a892a1SMatthew G. Knepley         PetscInt        minOrnt = -1;
212b5a892a1SMatthew G. Knepley         PetscInt        maxOrnt = 1;
213c4762a1bSJed Brown         PetscInt        **perms;
214c4762a1bSJed Brown 
2155f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscCalloc1(maxOrnt - minOrnt,&perms));
216c4762a1bSJed Brown         for (o = minOrnt; o < maxOrnt; o++) {
217c4762a1bSJed Brown           PetscInt *perm;
218c4762a1bSJed Brown 
219b5a892a1SMatthew G. Knepley           if (!o) { /* identity */
220c4762a1bSJed Brown             perms[o - minOrnt] = NULL;
221c4762a1bSJed Brown           } else {
2225f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscMalloc1(numDof * numComp, &perm));
223c4762a1bSJed Brown             for (i = numDof - 1, k = 0; i >= 0; i--) {
224c4762a1bSJed Brown               for (j = 0; j < numComp; j++, k++) perm[k] = i * numComp + j;
225c4762a1bSJed Brown             }
226c4762a1bSJed Brown             perms[o - minOrnt] = perm;
227c4762a1bSJed Brown           }
228c4762a1bSJed Brown         }
2295f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionSymLabelSetStratum(sym,d,numDof*numComp,minOrnt,maxOrnt,PETSC_OWN_POINTER,(const PetscInt **) perms,NULL));
230c4762a1bSJed Brown       } else if (d == 2) {
231c4762a1bSJed Brown         PetscInt        perEdge = user->k[f] - 1;
232c4762a1bSJed Brown         PetscInt        numDof  = perEdge * perEdge;
233c4762a1bSJed Brown         PetscInt        numComp = user->Nc[f];
234c4762a1bSJed Brown         PetscInt        minOrnt = -4;
235c4762a1bSJed Brown         PetscInt        maxOrnt = 4;
236c4762a1bSJed Brown         PetscInt        **perms;
237c4762a1bSJed Brown 
2385f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscCalloc1(maxOrnt-minOrnt,&perms));
239c4762a1bSJed Brown         for (o = minOrnt; o < maxOrnt; o++) {
240c4762a1bSJed Brown           PetscInt *perm;
241c4762a1bSJed Brown 
242c4762a1bSJed Brown           if (!o) continue; /* identity */
2435f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscMalloc1(numDof * numComp, &perm));
244c4762a1bSJed Brown           /* We want to perm[k] to list which *localArray* position the *sectionArray* position k should go to for the given orientation*/
245c4762a1bSJed Brown           switch (o) {
246c4762a1bSJed Brown           case 0:
247c4762a1bSJed Brown             break; /* identity */
248b5a892a1SMatthew G. Knepley           case -2: /* flip along (-1,-1)--( 1, 1), which swaps edges 0 and 3 and edges 1 and 2.  This swaps the i and j variables */
249c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
250c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
251c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
252c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * j + i) * numComp + c;
253c4762a1bSJed Brown                 }
254c4762a1bSJed Brown               }
255c4762a1bSJed Brown             }
256c4762a1bSJed Brown             break;
257b5a892a1SMatthew G. Knepley           case -1: /* flip along (-1, 0)--( 1, 0), which swaps edges 0 and 2.  This reverses the i variable */
258c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
259c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
260c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
261c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + j) * numComp + c;
262c4762a1bSJed Brown                 }
263c4762a1bSJed Brown               }
264c4762a1bSJed Brown             }
265c4762a1bSJed Brown             break;
266b5a892a1SMatthew G. Knepley           case -4: /* 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 */
267c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
268c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
269c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
270c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + (perEdge - 1 - i)) * numComp + c;
271c4762a1bSJed Brown                 }
272c4762a1bSJed Brown               }
273c4762a1bSJed Brown             }
274c4762a1bSJed Brown             break;
275b5a892a1SMatthew G. Knepley           case -3: /* flip along ( 0,-1)--( 0, 1), which swaps edges 3 and 1.  This reverses the j variable */
276c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
277c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
278c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
279c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * i + (perEdge - 1 - j)) * numComp + c;
280c4762a1bSJed Brown                 }
281c4762a1bSJed Brown               }
282c4762a1bSJed Brown             }
283c4762a1bSJed Brown             break;
284c4762a1bSJed Brown           case  1: /* rotate section edge 1 to local edge 0.  This swaps the i and j variables and then reverses the j variable */
285c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
286c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
287c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
288c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + i) * numComp + c;
289c4762a1bSJed Brown                 }
290c4762a1bSJed Brown               }
291c4762a1bSJed Brown             }
292c4762a1bSJed Brown             break;
293c4762a1bSJed Brown           case  2: /* rotate section edge 2 to local edge 0.  This reverse both i and j variables */
294c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
295c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
296c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
297c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + (perEdge - 1 - j)) * numComp + c;
298c4762a1bSJed Brown                 }
299c4762a1bSJed Brown               }
300c4762a1bSJed Brown             }
301c4762a1bSJed Brown             break;
302c4762a1bSJed Brown           case  3: /* rotate section edge 3 to local edge 0.  This swaps the i and j variables and then reverses the i variable */
303c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
304c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
305c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
306c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * j + (perEdge - 1 - i)) * numComp + c;
307c4762a1bSJed Brown                 }
308c4762a1bSJed Brown               }
309c4762a1bSJed Brown             }
310c4762a1bSJed Brown             break;
311c4762a1bSJed Brown           default:
312c4762a1bSJed Brown             break;
313c4762a1bSJed Brown           }
314c4762a1bSJed Brown           perms[o - minOrnt] = perm;
315c4762a1bSJed Brown         }
3165f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionSymLabelSetStratum(sym,d,numDof*numComp,minOrnt,maxOrnt,PETSC_OWN_POINTER,(const PetscInt **) perms,NULL));
317c4762a1bSJed Brown       }
318c4762a1bSJed Brown     }
3195f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetFieldSym(s,f,sym));
3205f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSymDestroy(&sym));
321c4762a1bSJed Brown   }
3225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionViewFromOptions(s,NULL,"-section_with_sym_view"));
323c4762a1bSJed Brown   PetscFunctionReturn(0);
324c4762a1bSJed Brown }
325c4762a1bSJed Brown 
326c4762a1bSJed Brown int main(int argc, char **argv)
327c4762a1bSJed Brown {
328c4762a1bSJed Brown   DM             dm;
329c4762a1bSJed Brown   PetscSection   s;
330c4762a1bSJed Brown   Vec            u;
331c4762a1bSJed Brown   AppCtx         user;
33230602db0SMatthew G. Knepley   PetscInt       dim, size = 0, f;
333c4762a1bSJed Brown 
334*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL,help));
3355f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
3365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm));
3375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(dm, DMPLEX));
3385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
3395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view"));
3405f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
341c4762a1bSJed Brown   /* Create a section for SEM order k */
342c4762a1bSJed Brown   {
343c4762a1bSJed Brown     PetscInt *numDof, d;
344c4762a1bSJed Brown 
3455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(user.Nf*(dim+1), &numDof));
346c4762a1bSJed Brown     for (f = 0; f < user.Nf; ++f) {
34730602db0SMatthew G. Knepley       for (d = 0; d <= dim; ++d) numDof[f*(dim+1)+d] = PetscPowInt(user.k[f]-1, d)*user.Nc[f];
348c4762a1bSJed Brown       size += PetscPowInt(user.k[f]+1, d)*user.Nc[f];
349c4762a1bSJed Brown     }
3505f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetNumFields(dm, user.Nf));
3515f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateSection(dm, NULL, user.Nc, numDof, 0, NULL, NULL, NULL, NULL, &s));
3525f80ce2aSJacob Faibussowitsch     CHKERRQ(SetSymmetries(dm, s, &user));
3535f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(numDof));
354c4762a1bSJed Brown   }
3555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetLocalSection(dm, s));
356c4762a1bSJed Brown   /* Create spectral ordering and load in data */
3575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
3585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(dm, &u));
35930602db0SMatthew G. Knepley   switch (dim) {
3605f80ce2aSJacob Faibussowitsch   case 2: CHKERRQ(LoadData2D(dm, 2, 2, size, u, &user));break;
3615f80ce2aSJacob Faibussowitsch   case 3: CHKERRQ(LoadData3D(dm, 2, 2, 2, size, u, &user));break;
362c4762a1bSJed Brown   }
363c4762a1bSJed Brown   /* Remove ordering and check some values */
3645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetClosurePermutation(s, (PetscObject) dm, dim, NULL));
36530602db0SMatthew G. Knepley   switch (dim) {
366c4762a1bSJed Brown   case 2:
3675f80ce2aSJacob Faibussowitsch     CHKERRQ(CheckPoint(dm, u,  0, &user));
3685f80ce2aSJacob Faibussowitsch     CHKERRQ(CheckPoint(dm, u, 13, &user));
3695f80ce2aSJacob Faibussowitsch     CHKERRQ(CheckPoint(dm, u, 15, &user));
3705f80ce2aSJacob Faibussowitsch     CHKERRQ(CheckPoint(dm, u, 19, &user));
371c4762a1bSJed Brown     break;
372c4762a1bSJed Brown   case 3:
3735f80ce2aSJacob Faibussowitsch     CHKERRQ(CheckPoint(dm, u,  0, &user));
3745f80ce2aSJacob Faibussowitsch     CHKERRQ(CheckPoint(dm, u, 13, &user));
3755f80ce2aSJacob Faibussowitsch     CHKERRQ(CheckPoint(dm, u, 15, &user));
3765f80ce2aSJacob Faibussowitsch     CHKERRQ(CheckPoint(dm, u, 19, &user));
377c4762a1bSJed Brown     break;
378c4762a1bSJed Brown   }
379c4762a1bSJed Brown   /* Recreate spectral ordering and read out data */
3805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s));
38130602db0SMatthew G. Knepley   switch (dim) {
3825f80ce2aSJacob Faibussowitsch   case 2: CHKERRQ(ReadData2D(dm, u, &user));break;
3835f80ce2aSJacob Faibussowitsch   case 3: CHKERRQ(ReadData3D(dm, u, &user));break;
384c4762a1bSJed Brown   }
3855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(dm, &u));
3865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&s));
3875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
3885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.Nc));
3895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.k));
390*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
391*b122ec5aSJacob Faibussowitsch   return 0;
392c4762a1bSJed Brown }
393c4762a1bSJed Brown 
394c4762a1bSJed Brown /*TEST
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   # Spectral ordering 2D 0-5
39730602db0SMatthew G. Knepley   testset:
39830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2
39930602db0SMatthew G. Knepley 
400c4762a1bSJed Brown     test:
401c4762a1bSJed Brown       suffix: 0
40230602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 2
403c4762a1bSJed Brown     test:
404c4762a1bSJed Brown       suffix: 1
40530602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 3
406c4762a1bSJed Brown     test:
407c4762a1bSJed Brown       suffix: 2
40830602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 5
409c4762a1bSJed Brown     test:
410c4762a1bSJed Brown       suffix: 3
41130602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 2 -order 2
412c4762a1bSJed Brown     test:
413c4762a1bSJed Brown       suffix: 4
41430602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,1 -order 2,2
415c4762a1bSJed Brown     test:
416c4762a1bSJed Brown       suffix: 5
41730602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,2 -order 2,3
41830602db0SMatthew G. Knepley 
419c4762a1bSJed Brown   # Spectral ordering 3D 6-11
42030602db0SMatthew G. Knepley   testset:
42130602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,2
42230602db0SMatthew G. Knepley 
423c4762a1bSJed Brown     test:
424c4762a1bSJed Brown       suffix: 6
42530602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 2
426c4762a1bSJed Brown     test:
427c4762a1bSJed Brown       suffix: 7
42830602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 3
429c4762a1bSJed Brown     test:
430c4762a1bSJed Brown       suffix: 8
43130602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 5
432c4762a1bSJed Brown     test:
433c4762a1bSJed Brown       suffix: 9
43430602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 2 -order 2
435c4762a1bSJed Brown     test:
436c4762a1bSJed Brown       suffix: 10
43730602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,1 -order 2,2
438c4762a1bSJed Brown     test:
439c4762a1bSJed Brown       suffix: 11
44030602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,2 -order 2,3
441c4762a1bSJed Brown 
442c4762a1bSJed Brown TEST*/
443