xref: /petsc/src/dm/impls/plex/tutorials/ex6.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
16c4762a1bSJed Brown   PetscFunctionBeginUser;
17c4762a1bSJed Brown   options->Nf = 0;
18c4762a1bSJed Brown   options->Nc = NULL;
19c4762a1bSJed Brown   options->k  = NULL;
20c4762a1bSJed Brown 
21*d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "SEM Problem Options", "DMPLEX");
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_fields", "The number of fields", "ex6.c", options->Nf, &options->Nf, NULL,0));
23c4762a1bSJed Brown   if (options->Nf) {
24c4762a1bSJed Brown     len  = options->Nf;
259566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, &options->Nc));
269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-num_components", "The number of components per field", "ex6.c", options->Nc, &len, &flg));
2708401ef6SPierre Jolivet     PetscCheck(!flg || !(len != options->Nf),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %d should be %d", len, options->Nf);
28c4762a1bSJed Brown     len  = options->Nf;
299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, &options->k));
309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-order", "The spectral order per field", "ex6.c", options->k, &len, &flg));
3108401ef6SPierre Jolivet     PetscCheck(!flg || !(len != options->Nf),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of order array is %d should be %d", len, options->Nf);
32c4762a1bSJed Brown   }
33*d0609cedSBarry Smith   PetscOptionsEnd();
34c4762a1bSJed Brown   PetscFunctionReturn(0);
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown static PetscErrorCode LoadData2D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt clSize, Vec u, AppCtx *user)
38c4762a1bSJed Brown {
39c4762a1bSJed Brown   PetscInt       i, j, f, c;
40c4762a1bSJed Brown   PetscScalar *closure;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   PetscFunctionBeginUser;
439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(clSize,&closure));
44c4762a1bSJed Brown   for (j = 0; j < Nj; ++j) {
45c4762a1bSJed Brown     for (i = 0; i < Ni; ++i) {
46c4762a1bSJed Brown       PetscInt    ki, kj, o = 0;
479566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(closure,clSize));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown       for (f = 0; f < user->Nf; ++f) {
50c4762a1bSJed Brown         PetscInt ioff = i*user->k[f], joff = j*user->k[f];
51c4762a1bSJed Brown 
52c4762a1bSJed Brown         for (kj = 0; kj <= user->k[f]; ++kj) {
53c4762a1bSJed Brown           for (ki = 0; ki <= user->k[f]; ++ki) {
54c4762a1bSJed Brown             for (c = 0; c < user->Nc[f]; ++c) {
55c4762a1bSJed Brown               closure[o++] = ((kj + joff)*(Ni*user->k[f]+1) + ki + ioff)*user->Nc[f]+c;
56c4762a1bSJed Brown             }
57c4762a1bSJed Brown           }
58c4762a1bSJed Brown         }
59c4762a1bSJed Brown       }
609566063dSJacob Faibussowitsch       PetscCall(DMPlexVecSetClosure(dm, NULL, u, j*Ni+i, closure, INSERT_VALUES));
61c4762a1bSJed Brown     }
62c4762a1bSJed Brown   }
639566063dSJacob Faibussowitsch   PetscCall(PetscFree(closure));
64c4762a1bSJed Brown   PetscFunctionReturn(0);
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown static PetscErrorCode LoadData3D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt Nk, PetscInt clSize, Vec u, AppCtx *user)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   PetscInt       i, j, k, f, c;
70c4762a1bSJed Brown   PetscScalar *closure;
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   PetscFunctionBeginUser;
739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(clSize,&closure));
74c4762a1bSJed Brown   for (k = 0; k < Nk; ++k) {
75c4762a1bSJed Brown     for (j = 0; j < Nj; ++j) {
76c4762a1bSJed Brown       for (i = 0; i < Ni; ++i) {
77c4762a1bSJed Brown         PetscInt    ki, kj, kk, o = 0;
789566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(closure,clSize));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown         for (f = 0; f < user->Nf; ++f) {
81c4762a1bSJed Brown           PetscInt ioff = i*user->k[f], joff = j*user->k[f], koff = k*user->k[f];
82c4762a1bSJed Brown 
83c4762a1bSJed Brown           for (kk = 0; kk <= user->k[f]; ++kk) {
84c4762a1bSJed Brown             for (kj = 0; kj <= user->k[f]; ++kj) {
85c4762a1bSJed Brown               for (ki = 0; ki <= user->k[f]; ++ki) {
86c4762a1bSJed Brown                 for (c = 0; c < user->Nc[f]; ++c) {
87c4762a1bSJed Brown                   closure[o++] = (((kk + koff)*(Nj*user->k[f]+1) + kj + joff)*(Ni*user->k[f]+1) + ki + ioff)*user->Nc[f]+c;
88c4762a1bSJed Brown                 }
89c4762a1bSJed Brown               }
90c4762a1bSJed Brown             }
91c4762a1bSJed Brown           }
92c4762a1bSJed Brown         }
939566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetClosure(dm, NULL, u, (k*Nj+j)*Ni+i, closure, INSERT_VALUES));
94c4762a1bSJed Brown       }
95c4762a1bSJed Brown     }
96c4762a1bSJed Brown   }
979566063dSJacob Faibussowitsch   PetscCall(PetscFree(closure));
98c4762a1bSJed Brown   PetscFunctionReturn(0);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown static PetscErrorCode CheckPoint(DM dm, Vec u, PetscInt point, AppCtx *user)
102c4762a1bSJed Brown {
103c4762a1bSJed Brown   PetscSection       s;
104c4762a1bSJed Brown   PetscScalar        *a;
105c4762a1bSJed Brown   const PetscScalar  *array;
106c4762a1bSJed Brown   PetscInt           dof, d;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   PetscFunctionBeginUser;
1099566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &s));
1109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &array));
1119566063dSJacob Faibussowitsch   PetscCall(DMPlexPointLocalRead(dm, point, array, &a));
1129566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetDof(s, point, &dof));
1139566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %D: ", point));
114c4762a1bSJed Brown   for (d = 0; d < dof; ++d) {
1159566063dSJacob Faibussowitsch     if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
1169566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double) PetscRealPart(a[d])));
117c4762a1bSJed Brown   }
1189566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
1199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &array));
120c4762a1bSJed Brown   PetscFunctionReturn(0);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown static PetscErrorCode ReadData2D(DM dm, Vec u, AppCtx *user)
124c4762a1bSJed Brown {
125c4762a1bSJed Brown   PetscInt       cStart, cEnd, cell;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   PetscFunctionBeginUser;
1289566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
129c4762a1bSJed Brown   for (cell = cStart; cell < cEnd; ++cell) {
130c4762a1bSJed Brown     PetscScalar *closure = NULL;
131c4762a1bSJed Brown     PetscInt     closureSize, ki, kj, f, c, foff = 0;
132c4762a1bSJed Brown 
1339566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure));
1349566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %D\n", cell));
135c4762a1bSJed Brown     for (f = 0; f < user->Nf; ++f) {
1369566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Field %D\n", f));
137c4762a1bSJed Brown       for (kj = user->k[f]; kj >= 0; --kj) {
138c4762a1bSJed Brown         for (ki = 0; ki <= user->k[f]; ++ki) {
1399566063dSJacob Faibussowitsch           if (ki > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  "));
140c4762a1bSJed Brown           for (c = 0; c < user->Nc[f]; ++c) {
1419566063dSJacob Faibussowitsch             if (c > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ","));
1429566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double) PetscRealPart(closure[(kj*(user->k[f]+1) + ki)*user->Nc[f]+c + foff])));
143c4762a1bSJed Brown           }
144c4762a1bSJed Brown         }
1459566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
146c4762a1bSJed Brown       }
1479566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
148c4762a1bSJed Brown       foff += PetscSqr(user->k[f]+1);
149c4762a1bSJed Brown     }
1509566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure));
1519566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
152c4762a1bSJed Brown   }
153c4762a1bSJed Brown   PetscFunctionReturn(0);
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown static PetscErrorCode ReadData3D(DM dm, Vec u, AppCtx *user)
157c4762a1bSJed Brown {
158c4762a1bSJed Brown   PetscInt       cStart, cEnd, cell;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   PetscFunctionBeginUser;
1619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
162c4762a1bSJed Brown   for (cell = cStart; cell < cEnd; ++cell) {
163c4762a1bSJed Brown     PetscScalar *closure = NULL;
164c4762a1bSJed Brown     PetscInt     closureSize, ki, kj, kk, f, c, foff = 0;
165c4762a1bSJed Brown 
1669566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure));
1679566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %D\n", cell));
168c4762a1bSJed Brown     for (f = 0; f < user->Nf; ++f) {
1699566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Field %D\n", f));
170c4762a1bSJed Brown       for (kk = user->k[f]; kk >= 0; --kk) {
171c4762a1bSJed Brown         for (kj = user->k[f]; kj >= 0; --kj) {
172c4762a1bSJed Brown           for (ki = 0; ki <= user->k[f]; ++ki) {
1739566063dSJacob Faibussowitsch             if (ki > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  "));
174c4762a1bSJed Brown             for (c = 0; c < user->Nc[f]; ++c) {
1759566063dSJacob Faibussowitsch               if (c > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ","));
1769566063dSJacob Faibussowitsch               PetscCall(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])));
177c4762a1bSJed Brown             }
178c4762a1bSJed Brown           }
1799566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
180c4762a1bSJed Brown         }
1819566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
182c4762a1bSJed Brown       }
1839566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
184c4762a1bSJed Brown       foff += PetscSqr(user->k[f]+1);
185c4762a1bSJed Brown     }
1869566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure));
1879566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
188c4762a1bSJed Brown   }
189c4762a1bSJed Brown   PetscFunctionReturn(0);
190c4762a1bSJed Brown }
191c4762a1bSJed Brown 
192c4762a1bSJed Brown static PetscErrorCode SetSymmetries(DM dm, PetscSection s, AppCtx *user)
193c4762a1bSJed Brown {
19430602db0SMatthew G. Knepley   PetscInt       dim, f, o, i, j, k, c, d;
195c4762a1bSJed Brown   DMLabel        depthLabel;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   PetscFunctionBegin;
1989566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1999566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm,"depth",&depthLabel));
200c4762a1bSJed Brown   for (f = 0; f < user->Nf; f++) {
201c4762a1bSJed Brown     PetscSectionSym sym;
202c4762a1bSJed Brown 
203c4762a1bSJed 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 */
2049566063dSJacob Faibussowitsch     PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)s),depthLabel,&sym));
205c4762a1bSJed Brown 
20630602db0SMatthew G. Knepley     for (d = 0; d <= dim; d++) {
207c4762a1bSJed Brown       if (d == 1) {
208c4762a1bSJed Brown         PetscInt        numDof  = user->k[f] - 1;
209c4762a1bSJed Brown         PetscInt        numComp = user->Nc[f];
210b5a892a1SMatthew G. Knepley         PetscInt        minOrnt = -1;
211b5a892a1SMatthew G. Knepley         PetscInt        maxOrnt = 1;
212c4762a1bSJed Brown         PetscInt        **perms;
213c4762a1bSJed Brown 
2149566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(maxOrnt - minOrnt,&perms));
215c4762a1bSJed Brown         for (o = minOrnt; o < maxOrnt; o++) {
216c4762a1bSJed Brown           PetscInt *perm;
217c4762a1bSJed Brown 
218b5a892a1SMatthew G. Knepley           if (!o) { /* identity */
219c4762a1bSJed Brown             perms[o - minOrnt] = NULL;
220c4762a1bSJed Brown           } else {
2219566063dSJacob Faibussowitsch             PetscCall(PetscMalloc1(numDof * numComp, &perm));
222c4762a1bSJed Brown             for (i = numDof - 1, k = 0; i >= 0; i--) {
223c4762a1bSJed Brown               for (j = 0; j < numComp; j++, k++) perm[k] = i * numComp + j;
224c4762a1bSJed Brown             }
225c4762a1bSJed Brown             perms[o - minOrnt] = perm;
226c4762a1bSJed Brown           }
227c4762a1bSJed Brown         }
2289566063dSJacob Faibussowitsch         PetscCall(PetscSectionSymLabelSetStratum(sym,d,numDof*numComp,minOrnt,maxOrnt,PETSC_OWN_POINTER,(const PetscInt **) perms,NULL));
229c4762a1bSJed Brown       } else if (d == 2) {
230c4762a1bSJed Brown         PetscInt        perEdge = user->k[f] - 1;
231c4762a1bSJed Brown         PetscInt        numDof  = perEdge * perEdge;
232c4762a1bSJed Brown         PetscInt        numComp = user->Nc[f];
233c4762a1bSJed Brown         PetscInt        minOrnt = -4;
234c4762a1bSJed Brown         PetscInt        maxOrnt = 4;
235c4762a1bSJed Brown         PetscInt        **perms;
236c4762a1bSJed Brown 
2379566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(maxOrnt-minOrnt,&perms));
238c4762a1bSJed Brown         for (o = minOrnt; o < maxOrnt; o++) {
239c4762a1bSJed Brown           PetscInt *perm;
240c4762a1bSJed Brown 
241c4762a1bSJed Brown           if (!o) continue; /* identity */
2429566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(numDof * numComp, &perm));
243c4762a1bSJed Brown           /* We want to perm[k] to list which *localArray* position the *sectionArray* position k should go to for the given orientation*/
244c4762a1bSJed Brown           switch (o) {
245c4762a1bSJed Brown           case 0:
246c4762a1bSJed Brown             break; /* identity */
247b5a892a1SMatthew 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 */
248c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
249c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
250c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
251c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * j + i) * numComp + c;
252c4762a1bSJed Brown                 }
253c4762a1bSJed Brown               }
254c4762a1bSJed Brown             }
255c4762a1bSJed Brown             break;
256b5a892a1SMatthew G. Knepley           case -1: /* flip along (-1, 0)--( 1, 0), which swaps edges 0 and 2.  This reverses the i variable */
257c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
258c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
259c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
260c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + j) * numComp + c;
261c4762a1bSJed Brown                 }
262c4762a1bSJed Brown               }
263c4762a1bSJed Brown             }
264c4762a1bSJed Brown             break;
265b5a892a1SMatthew 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 */
266c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
267c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
268c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
269c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + (perEdge - 1 - i)) * numComp + c;
270c4762a1bSJed Brown                 }
271c4762a1bSJed Brown               }
272c4762a1bSJed Brown             }
273c4762a1bSJed Brown             break;
274b5a892a1SMatthew G. Knepley           case -3: /* flip along ( 0,-1)--( 0, 1), which swaps edges 3 and 1.  This reverses the j variable */
275c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
276c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
277c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
278c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * i + (perEdge - 1 - j)) * numComp + c;
279c4762a1bSJed Brown                 }
280c4762a1bSJed Brown               }
281c4762a1bSJed Brown             }
282c4762a1bSJed Brown             break;
283c4762a1bSJed Brown           case  1: /* rotate section edge 1 to local edge 0.  This swaps the i and j variables and then reverses the j variable */
284c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
285c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
286c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
287c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + i) * numComp + c;
288c4762a1bSJed Brown                 }
289c4762a1bSJed Brown               }
290c4762a1bSJed Brown             }
291c4762a1bSJed Brown             break;
292c4762a1bSJed Brown           case  2: /* rotate section edge 2 to local edge 0.  This reverse both i and j variables */
293c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
294c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
295c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
296c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + (perEdge - 1 - j)) * numComp + c;
297c4762a1bSJed Brown                 }
298c4762a1bSJed Brown               }
299c4762a1bSJed Brown             }
300c4762a1bSJed Brown             break;
301c4762a1bSJed Brown           case  3: /* rotate section edge 3 to local edge 0.  This swaps the i and j variables and then reverses the i variable */
302c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
303c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
304c4762a1bSJed Brown                 for (c = 0; c < numComp; c++) {
305c4762a1bSJed Brown                   perm[k * numComp + c] = (perEdge * j + (perEdge - 1 - i)) * numComp + c;
306c4762a1bSJed Brown                 }
307c4762a1bSJed Brown               }
308c4762a1bSJed Brown             }
309c4762a1bSJed Brown             break;
310c4762a1bSJed Brown           default:
311c4762a1bSJed Brown             break;
312c4762a1bSJed Brown           }
313c4762a1bSJed Brown           perms[o - minOrnt] = perm;
314c4762a1bSJed Brown         }
3159566063dSJacob Faibussowitsch         PetscCall(PetscSectionSymLabelSetStratum(sym,d,numDof*numComp,minOrnt,maxOrnt,PETSC_OWN_POINTER,(const PetscInt **) perms,NULL));
316c4762a1bSJed Brown       }
317c4762a1bSJed Brown     }
3189566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldSym(s,f,sym));
3199566063dSJacob Faibussowitsch     PetscCall(PetscSectionSymDestroy(&sym));
320c4762a1bSJed Brown   }
3219566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(s,NULL,"-section_with_sym_view"));
322c4762a1bSJed Brown   PetscFunctionReturn(0);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown 
325c4762a1bSJed Brown int main(int argc, char **argv)
326c4762a1bSJed Brown {
327c4762a1bSJed Brown   DM             dm;
328c4762a1bSJed Brown   PetscSection   s;
329c4762a1bSJed Brown   Vec            u;
330c4762a1bSJed Brown   AppCtx         user;
33130602db0SMatthew G. Knepley   PetscInt       dim, size = 0, f;
332c4762a1bSJed Brown 
3339566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
3349566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3359566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
3369566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
3379566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
3389566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
3399566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
340c4762a1bSJed Brown   /* Create a section for SEM order k */
341c4762a1bSJed Brown   {
342c4762a1bSJed Brown     PetscInt *numDof, d;
343c4762a1bSJed Brown 
3449566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(user.Nf*(dim+1), &numDof));
345c4762a1bSJed Brown     for (f = 0; f < user.Nf; ++f) {
34630602db0SMatthew G. Knepley       for (d = 0; d <= dim; ++d) numDof[f*(dim+1)+d] = PetscPowInt(user.k[f]-1, d)*user.Nc[f];
347c4762a1bSJed Brown       size += PetscPowInt(user.k[f]+1, d)*user.Nc[f];
348c4762a1bSJed Brown     }
3499566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(dm, user.Nf));
3509566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(dm, NULL, user.Nc, numDof, 0, NULL, NULL, NULL, NULL, &s));
3519566063dSJacob Faibussowitsch     PetscCall(SetSymmetries(dm, s, &user));
3529566063dSJacob Faibussowitsch     PetscCall(PetscFree(numDof));
353c4762a1bSJed Brown   }
3549566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, s));
355c4762a1bSJed Brown   /* Create spectral ordering and load in data */
3569566063dSJacob Faibussowitsch   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
3579566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &u));
35830602db0SMatthew G. Knepley   switch (dim) {
3599566063dSJacob Faibussowitsch   case 2: PetscCall(LoadData2D(dm, 2, 2, size, u, &user));break;
3609566063dSJacob Faibussowitsch   case 3: PetscCall(LoadData3D(dm, 2, 2, 2, size, u, &user));break;
361c4762a1bSJed Brown   }
362c4762a1bSJed Brown   /* Remove ordering and check some values */
3639566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetClosurePermutation(s, (PetscObject) dm, dim, NULL));
36430602db0SMatthew G. Knepley   switch (dim) {
365c4762a1bSJed Brown   case 2:
3669566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u,  0, &user));
3679566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 13, &user));
3689566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 15, &user));
3699566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 19, &user));
370c4762a1bSJed Brown     break;
371c4762a1bSJed Brown   case 3:
3729566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u,  0, &user));
3739566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 13, &user));
3749566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 15, &user));
3759566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 19, &user));
376c4762a1bSJed Brown     break;
377c4762a1bSJed Brown   }
378c4762a1bSJed Brown   /* Recreate spectral ordering and read out data */
3799566063dSJacob Faibussowitsch   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s));
38030602db0SMatthew G. Knepley   switch (dim) {
3819566063dSJacob Faibussowitsch   case 2: PetscCall(ReadData2D(dm, u, &user));break;
3829566063dSJacob Faibussowitsch   case 3: PetscCall(ReadData3D(dm, u, &user));break;
383c4762a1bSJed Brown   }
3849566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &u));
3859566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s));
3869566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3879566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.Nc));
3889566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.k));
3899566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
390b122ec5aSJacob Faibussowitsch   return 0;
391c4762a1bSJed Brown }
392c4762a1bSJed Brown 
393c4762a1bSJed Brown /*TEST
394c4762a1bSJed Brown 
395c4762a1bSJed Brown   # Spectral ordering 2D 0-5
39630602db0SMatthew G. Knepley   testset:
39730602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2
39830602db0SMatthew G. Knepley 
399c4762a1bSJed Brown     test:
400c4762a1bSJed Brown       suffix: 0
40130602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 2
402c4762a1bSJed Brown     test:
403c4762a1bSJed Brown       suffix: 1
40430602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 3
405c4762a1bSJed Brown     test:
406c4762a1bSJed Brown       suffix: 2
40730602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 5
408c4762a1bSJed Brown     test:
409c4762a1bSJed Brown       suffix: 3
41030602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 2 -order 2
411c4762a1bSJed Brown     test:
412c4762a1bSJed Brown       suffix: 4
41330602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,1 -order 2,2
414c4762a1bSJed Brown     test:
415c4762a1bSJed Brown       suffix: 5
41630602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,2 -order 2,3
41730602db0SMatthew G. Knepley 
418c4762a1bSJed Brown   # Spectral ordering 3D 6-11
41930602db0SMatthew G. Knepley   testset:
42030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,2
42130602db0SMatthew G. Knepley 
422c4762a1bSJed Brown     test:
423c4762a1bSJed Brown       suffix: 6
42430602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 2
425c4762a1bSJed Brown     test:
426c4762a1bSJed Brown       suffix: 7
42730602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 3
428c4762a1bSJed Brown     test:
429c4762a1bSJed Brown       suffix: 8
43030602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 5
431c4762a1bSJed Brown     test:
432c4762a1bSJed Brown       suffix: 9
43330602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 2 -order 2
434c4762a1bSJed Brown     test:
435c4762a1bSJed Brown       suffix: 10
43630602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,1 -order 2,2
437c4762a1bSJed Brown     test:
438c4762a1bSJed Brown       suffix: 11
43930602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,2 -order 2,3
440c4762a1bSJed Brown 
441c4762a1bSJed Brown TEST*/
442