xref: /petsc/src/dm/impls/plex/tutorials/ex6.c (revision 31d78bcd2b98084dc1368b20eb1129c8b9fb39fe)
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 
ProcessOptions(MPI_Comm comm,AppCtx * options)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12d71ae5a4SJacob Faibussowitsch {
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 
21d0609cedSBarry 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));
2763a3b9bcSJacob Faibussowitsch     PetscCheck(!flg || !(len != options->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %" PetscInt_FMT " should be %" PetscInt_FMT, 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));
3163a3b9bcSJacob Faibussowitsch     PetscCheck(!flg || !(len != options->Nf), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of order array is %" PetscInt_FMT " should be %" PetscInt_FMT, len, options->Nf);
32c4762a1bSJed Brown   }
33d0609cedSBarry Smith   PetscOptionsEnd();
34*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
LoadData2D(DM dm,PetscInt Ni,PetscInt Nj,PetscInt clSize,Vec u,AppCtx * user)37d71ae5a4SJacob Faibussowitsch static PetscErrorCode LoadData2D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt clSize, Vec u, AppCtx *user)
38d71ae5a4SJacob Faibussowitsch {
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) {
54ad540459SPierre Jolivet             for (c = 0; c < user->Nc[f]; ++c) closure[o++] = ((kj + joff) * (Ni * user->k[f] + 1) + ki + ioff) * user->Nc[f] + c;
55c4762a1bSJed Brown           }
56c4762a1bSJed Brown         }
57c4762a1bSJed Brown       }
589566063dSJacob Faibussowitsch       PetscCall(DMPlexVecSetClosure(dm, NULL, u, j * Ni + i, closure, INSERT_VALUES));
59c4762a1bSJed Brown     }
60c4762a1bSJed Brown   }
619566063dSJacob Faibussowitsch   PetscCall(PetscFree(closure));
62*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63c4762a1bSJed Brown }
64c4762a1bSJed Brown 
LoadData3D(DM dm,PetscInt Ni,PetscInt Nj,PetscInt Nk,PetscInt clSize,Vec u,AppCtx * user)65d71ae5a4SJacob Faibussowitsch static PetscErrorCode LoadData3D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt Nk, PetscInt clSize, Vec u, AppCtx *user)
66d71ae5a4SJacob Faibussowitsch {
67c4762a1bSJed Brown   PetscInt     i, j, k, f, c;
68c4762a1bSJed Brown   PetscScalar *closure;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBeginUser;
719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(clSize, &closure));
72c4762a1bSJed Brown   for (k = 0; k < Nk; ++k) {
73c4762a1bSJed Brown     for (j = 0; j < Nj; ++j) {
74c4762a1bSJed Brown       for (i = 0; i < Ni; ++i) {
75c4762a1bSJed Brown         PetscInt ki, kj, kk, o = 0;
769566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(closure, clSize));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown         for (f = 0; f < user->Nf; ++f) {
79c4762a1bSJed Brown           PetscInt ioff = i * user->k[f], joff = j * user->k[f], koff = k * user->k[f];
80c4762a1bSJed Brown 
81c4762a1bSJed Brown           for (kk = 0; kk <= user->k[f]; ++kk) {
82c4762a1bSJed Brown             for (kj = 0; kj <= user->k[f]; ++kj) {
83c4762a1bSJed Brown               for (ki = 0; ki <= user->k[f]; ++ki) {
84ad540459SPierre Jolivet                 for (c = 0; c < user->Nc[f]; ++c) closure[o++] = (((kk + koff) * (Nj * user->k[f] + 1) + kj + joff) * (Ni * user->k[f] + 1) + ki + ioff) * user->Nc[f] + c;
85c4762a1bSJed Brown               }
86c4762a1bSJed Brown             }
87c4762a1bSJed Brown           }
88c4762a1bSJed Brown         }
899566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetClosure(dm, NULL, u, (k * Nj + j) * Ni + i, closure, INSERT_VALUES));
90c4762a1bSJed Brown       }
91c4762a1bSJed Brown     }
92c4762a1bSJed Brown   }
939566063dSJacob Faibussowitsch   PetscCall(PetscFree(closure));
94*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
CheckPoint(DM dm,Vec u,PetscInt point,AppCtx * user)97d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckPoint(DM dm, Vec u, PetscInt point, AppCtx *user)
98d71ae5a4SJacob Faibussowitsch {
99c4762a1bSJed Brown   PetscSection       s;
100c4762a1bSJed Brown   PetscScalar       *a;
101c4762a1bSJed Brown   const PetscScalar *array;
102c4762a1bSJed Brown   PetscInt           dof, d;
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   PetscFunctionBeginUser;
1059566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &s));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(u, &array));
1079566063dSJacob Faibussowitsch   PetscCall(DMPlexPointLocalRead(dm, point, array, &a));
1089566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetDof(s, point, &dof));
10963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT ": ", point));
110c4762a1bSJed Brown   for (d = 0; d < dof; ++d) {
1119566063dSJacob Faibussowitsch     if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
1129566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double)PetscRealPart(a[d])));
113c4762a1bSJed Brown   }
1149566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
1159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(u, &array));
116*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117c4762a1bSJed Brown }
118c4762a1bSJed Brown 
ReadData2D(DM dm,Vec u,AppCtx * user)119d71ae5a4SJacob Faibussowitsch static PetscErrorCode ReadData2D(DM dm, Vec u, AppCtx *user)
120d71ae5a4SJacob Faibussowitsch {
121c4762a1bSJed Brown   PetscInt cStart, cEnd, cell;
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   PetscFunctionBeginUser;
1249566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
125c4762a1bSJed Brown   for (cell = cStart; cell < cEnd; ++cell) {
126c4762a1bSJed Brown     PetscScalar *closure = NULL;
127c4762a1bSJed Brown     PetscInt     closureSize, ki, kj, f, c, foff = 0;
128c4762a1bSJed Brown 
1299566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure));
13063a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT "\n", cell));
131c4762a1bSJed Brown     for (f = 0; f < user->Nf; ++f) {
13263a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Field %" PetscInt_FMT "\n", f));
133c4762a1bSJed Brown       for (kj = user->k[f]; kj >= 0; --kj) {
134c4762a1bSJed Brown         for (ki = 0; ki <= user->k[f]; ++ki) {
1359566063dSJacob Faibussowitsch           if (ki > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  "));
136c4762a1bSJed Brown           for (c = 0; c < user->Nc[f]; ++c) {
1379566063dSJacob Faibussowitsch             if (c > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ","));
1389566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double)PetscRealPart(closure[(kj * (user->k[f] + 1) + ki) * user->Nc[f] + c + foff])));
139c4762a1bSJed Brown           }
140c4762a1bSJed Brown         }
1419566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
142c4762a1bSJed Brown       }
1439566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
144c4762a1bSJed Brown       foff += PetscSqr(user->k[f] + 1);
145c4762a1bSJed Brown     }
1469566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure));
1479566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
148c4762a1bSJed Brown   }
149*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
ReadData3D(DM dm,Vec u,AppCtx * user)152d71ae5a4SJacob Faibussowitsch static PetscErrorCode ReadData3D(DM dm, Vec u, AppCtx *user)
153d71ae5a4SJacob Faibussowitsch {
154c4762a1bSJed Brown   PetscInt cStart, cEnd, cell;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   PetscFunctionBeginUser;
1579566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
158c4762a1bSJed Brown   for (cell = cStart; cell < cEnd; ++cell) {
159c4762a1bSJed Brown     PetscScalar *closure = NULL;
160c4762a1bSJed Brown     PetscInt     closureSize, ki, kj, kk, f, c, foff = 0;
161c4762a1bSJed Brown 
1629566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure));
16363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT "\n", cell));
164c4762a1bSJed Brown     for (f = 0; f < user->Nf; ++f) {
16563a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Field %" PetscInt_FMT "\n", f));
166c4762a1bSJed Brown       for (kk = user->k[f]; kk >= 0; --kk) {
167c4762a1bSJed Brown         for (kj = user->k[f]; kj >= 0; --kj) {
168c4762a1bSJed Brown           for (ki = 0; ki <= user->k[f]; ++ki) {
1699566063dSJacob Faibussowitsch             if (ki > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  "));
170c4762a1bSJed Brown             for (c = 0; c < user->Nc[f]; ++c) {
1719566063dSJacob Faibussowitsch               if (c > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ","));
1729566063dSJacob 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])));
173c4762a1bSJed Brown             }
174c4762a1bSJed Brown           }
1759566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
176c4762a1bSJed Brown         }
1779566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
178c4762a1bSJed Brown       }
1799566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
180c4762a1bSJed Brown       foff += PetscSqr(user->k[f] + 1);
181c4762a1bSJed Brown     }
1829566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure));
1839566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n"));
184c4762a1bSJed Brown   }
185*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
SetSymmetries(DM dm,PetscSection s,AppCtx * user)188d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetSymmetries(DM dm, PetscSection s, AppCtx *user)
189d71ae5a4SJacob Faibussowitsch {
19030602db0SMatthew G. Knepley   PetscInt dim, f, o, i, j, k, c, d;
191c4762a1bSJed Brown   DMLabel  depthLabel;
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   PetscFunctionBegin;
1949566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1959566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "depth", &depthLabel));
196c4762a1bSJed Brown   for (f = 0; f < user->Nf; f++) {
197c4762a1bSJed Brown     PetscSectionSym sym;
198c4762a1bSJed Brown 
199c4762a1bSJed 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 */
2009566063dSJacob Faibussowitsch     PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)s), depthLabel, &sym));
201c4762a1bSJed Brown 
20230602db0SMatthew G. Knepley     for (d = 0; d <= dim; d++) {
203c4762a1bSJed Brown       if (d == 1) {
204c4762a1bSJed Brown         PetscInt   numDof  = user->k[f] - 1;
205c4762a1bSJed Brown         PetscInt   numComp = user->Nc[f];
206b5a892a1SMatthew G. Knepley         PetscInt   minOrnt = -1;
207b5a892a1SMatthew G. Knepley         PetscInt   maxOrnt = 1;
208c4762a1bSJed Brown         PetscInt **perms;
209c4762a1bSJed Brown 
2109566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(maxOrnt - minOrnt, &perms));
211c4762a1bSJed Brown         for (o = minOrnt; o < maxOrnt; o++) {
212c4762a1bSJed Brown           PetscInt *perm;
213c4762a1bSJed Brown 
214b5a892a1SMatthew G. Knepley           if (!o) { /* identity */
215c4762a1bSJed Brown             perms[o - minOrnt] = NULL;
216c4762a1bSJed Brown           } else {
2179566063dSJacob Faibussowitsch             PetscCall(PetscMalloc1(numDof * numComp, &perm));
218c4762a1bSJed Brown             for (i = numDof - 1, k = 0; i >= 0; i--) {
219c4762a1bSJed Brown               for (j = 0; j < numComp; j++, k++) perm[k] = i * numComp + j;
220c4762a1bSJed Brown             }
221c4762a1bSJed Brown             perms[o - minOrnt] = perm;
222c4762a1bSJed Brown           }
223c4762a1bSJed Brown         }
2249566063dSJacob Faibussowitsch         PetscCall(PetscSectionSymLabelSetStratum(sym, d, numDof * numComp, minOrnt, maxOrnt, PETSC_OWN_POINTER, (const PetscInt **)perms, NULL));
225c4762a1bSJed Brown       } else if (d == 2) {
226c4762a1bSJed Brown         PetscInt   perEdge = user->k[f] - 1;
227c4762a1bSJed Brown         PetscInt   numDof  = perEdge * perEdge;
228c4762a1bSJed Brown         PetscInt   numComp = user->Nc[f];
229c4762a1bSJed Brown         PetscInt   minOrnt = -4;
230c4762a1bSJed Brown         PetscInt   maxOrnt = 4;
231c4762a1bSJed Brown         PetscInt **perms;
232c4762a1bSJed Brown 
2339566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(maxOrnt - minOrnt, &perms));
234c4762a1bSJed Brown         for (o = minOrnt; o < maxOrnt; o++) {
235c4762a1bSJed Brown           PetscInt *perm;
236c4762a1bSJed Brown 
237c4762a1bSJed Brown           if (!o) continue; /* identity */
2389566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(numDof * numComp, &perm));
239c4762a1bSJed Brown           /* We want to perm[k] to list which *localArray* position the *sectionArray* position k should go to for the given orientation*/
240c4762a1bSJed Brown           switch (o) {
241d71ae5a4SJacob Faibussowitsch           case 0:
242d71ae5a4SJacob Faibussowitsch             break; /* identity */
243b5a892a1SMatthew 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 */
244c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
245c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
246ad540459SPierre Jolivet                 for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * j + i) * numComp + c;
247c4762a1bSJed Brown               }
248c4762a1bSJed Brown             }
249c4762a1bSJed Brown             break;
250b5a892a1SMatthew G. Knepley           case -1: /* flip along (-1, 0)--( 1, 0), which swaps edges 0 and 2.  This reverses the i variable */
251c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
252c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
253ad540459SPierre Jolivet                 for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + j) * numComp + c;
254c4762a1bSJed Brown               }
255c4762a1bSJed Brown             }
256c4762a1bSJed Brown             break;
257b5a892a1SMatthew 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 */
258c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
259c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
260ad540459SPierre Jolivet                 for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + (perEdge - 1 - i)) * numComp + c;
261c4762a1bSJed Brown               }
262c4762a1bSJed Brown             }
263c4762a1bSJed Brown             break;
264b5a892a1SMatthew G. Knepley           case -3: /* flip along ( 0,-1)--( 0, 1), which swaps edges 3 and 1.  This reverses the j variable */
265c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
266c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
267ad540459SPierre Jolivet                 for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * i + (perEdge - 1 - j)) * numComp + c;
268c4762a1bSJed Brown               }
269c4762a1bSJed Brown             }
270c4762a1bSJed Brown             break;
271c4762a1bSJed Brown           case 1: /* rotate section edge 1 to local edge 0.  This swaps the i and j variables and then reverses the j variable */
272c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
273c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
274ad540459SPierre Jolivet                 for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + i) * numComp + c;
275c4762a1bSJed Brown               }
276c4762a1bSJed Brown             }
277c4762a1bSJed Brown             break;
278c4762a1bSJed Brown           case 2: /* rotate section edge 2 to local edge 0.  This reverse both i and j variables */
279c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
280c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
281ad540459SPierre Jolivet                 for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + (perEdge - 1 - j)) * numComp + c;
282c4762a1bSJed Brown               }
283c4762a1bSJed Brown             }
284c4762a1bSJed Brown             break;
285c4762a1bSJed Brown           case 3: /* rotate section edge 3 to local edge 0.  This swaps the i and j variables and then reverses the i variable */
286c4762a1bSJed Brown             for (i = 0, k = 0; i < perEdge; i++) {
287c4762a1bSJed Brown               for (j = 0; j < perEdge; j++, k++) {
288ad540459SPierre Jolivet                 for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * j + (perEdge - 1 - i)) * numComp + c;
289c4762a1bSJed Brown               }
290c4762a1bSJed Brown             }
291c4762a1bSJed Brown             break;
292d71ae5a4SJacob Faibussowitsch           default:
293d71ae5a4SJacob Faibussowitsch             break;
294c4762a1bSJed Brown           }
295c4762a1bSJed Brown           perms[o - minOrnt] = perm;
296c4762a1bSJed Brown         }
2979566063dSJacob Faibussowitsch         PetscCall(PetscSectionSymLabelSetStratum(sym, d, numDof * numComp, minOrnt, maxOrnt, PETSC_OWN_POINTER, (const PetscInt **)perms, NULL));
298c4762a1bSJed Brown       }
299c4762a1bSJed Brown     }
3009566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldSym(s, f, sym));
3019566063dSJacob Faibussowitsch     PetscCall(PetscSectionSymDestroy(&sym));
302c4762a1bSJed Brown   }
3039566063dSJacob Faibussowitsch   PetscCall(PetscSectionViewFromOptions(s, NULL, "-section_with_sym_view"));
304*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
305c4762a1bSJed Brown }
306c4762a1bSJed Brown 
main(int argc,char ** argv)307d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
308d71ae5a4SJacob Faibussowitsch {
309c4762a1bSJed Brown   DM           dm;
310c4762a1bSJed Brown   PetscSection s;
311c4762a1bSJed Brown   Vec          u;
312c4762a1bSJed Brown   AppCtx       user;
31330602db0SMatthew G. Knepley   PetscInt     dim, size = 0, f;
314c4762a1bSJed Brown 
315327415f7SBarry Smith   PetscFunctionBeginUser;
3169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3179566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3189566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
3199566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
3209566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
3219566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
3229566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
323c4762a1bSJed Brown   /* Create a section for SEM order k */
324c4762a1bSJed Brown   {
325c4762a1bSJed Brown     PetscInt *numDof, d;
326c4762a1bSJed Brown 
3279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(user.Nf * (dim + 1), &numDof));
328c4762a1bSJed Brown     for (f = 0; f < user.Nf; ++f) {
32930602db0SMatthew G. Knepley       for (d = 0; d <= dim; ++d) numDof[f * (dim + 1) + d] = PetscPowInt(user.k[f] - 1, d) * user.Nc[f];
330c4762a1bSJed Brown       size += PetscPowInt(user.k[f] + 1, d) * user.Nc[f];
331c4762a1bSJed Brown     }
3329566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(dm, user.Nf));
3339566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(dm, NULL, user.Nc, numDof, 0, NULL, NULL, NULL, NULL, &s));
3349566063dSJacob Faibussowitsch     PetscCall(SetSymmetries(dm, s, &user));
3359566063dSJacob Faibussowitsch     PetscCall(PetscFree(numDof));
336c4762a1bSJed Brown   }
3379566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, s));
338c4762a1bSJed Brown   /* Create spectral ordering and load in data */
3399566063dSJacob Faibussowitsch   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL));
3409566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &u));
34130602db0SMatthew G. Knepley   switch (dim) {
342d71ae5a4SJacob Faibussowitsch   case 2:
343d71ae5a4SJacob Faibussowitsch     PetscCall(LoadData2D(dm, 2, 2, size, u, &user));
344d71ae5a4SJacob Faibussowitsch     break;
345d71ae5a4SJacob Faibussowitsch   case 3:
346d71ae5a4SJacob Faibussowitsch     PetscCall(LoadData3D(dm, 2, 2, 2, size, u, &user));
347d71ae5a4SJacob Faibussowitsch     break;
348c4762a1bSJed Brown   }
349c4762a1bSJed Brown   /* Remove ordering and check some values */
3509566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetClosurePermutation(s, (PetscObject)dm, dim, NULL));
35130602db0SMatthew G. Knepley   switch (dim) {
352c4762a1bSJed Brown   case 2:
3539566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 0, &user));
3549566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 13, &user));
3559566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 15, &user));
3569566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 19, &user));
357c4762a1bSJed Brown     break;
358c4762a1bSJed Brown   case 3:
3599566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 0, &user));
3609566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 13, &user));
3619566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 15, &user));
3629566063dSJacob Faibussowitsch     PetscCall(CheckPoint(dm, u, 19, &user));
363c4762a1bSJed Brown     break;
364c4762a1bSJed Brown   }
365c4762a1bSJed Brown   /* Recreate spectral ordering and read out data */
3669566063dSJacob Faibussowitsch   PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s));
36730602db0SMatthew G. Knepley   switch (dim) {
368d71ae5a4SJacob Faibussowitsch   case 2:
369d71ae5a4SJacob Faibussowitsch     PetscCall(ReadData2D(dm, u, &user));
370d71ae5a4SJacob Faibussowitsch     break;
371d71ae5a4SJacob Faibussowitsch   case 3:
372d71ae5a4SJacob Faibussowitsch     PetscCall(ReadData3D(dm, u, &user));
373d71ae5a4SJacob Faibussowitsch     break;
374c4762a1bSJed Brown   }
3759566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &u));
3769566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&s));
3779566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3789566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.Nc));
3799566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.k));
3809566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
381b122ec5aSJacob Faibussowitsch   return 0;
382c4762a1bSJed Brown }
383c4762a1bSJed Brown 
384c4762a1bSJed Brown /*TEST
385c4762a1bSJed Brown 
386c4762a1bSJed Brown   # Spectral ordering 2D 0-5
38730602db0SMatthew G. Knepley   testset:
38830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2
38930602db0SMatthew G. Knepley 
390c4762a1bSJed Brown     test:
391c4762a1bSJed Brown       suffix: 0
39230602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 2
393c4762a1bSJed Brown     test:
394c4762a1bSJed Brown       suffix: 1
39530602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 3
396c4762a1bSJed Brown     test:
397c4762a1bSJed Brown       suffix: 2
39830602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 5
399c4762a1bSJed Brown     test:
400c4762a1bSJed Brown       suffix: 3
40130602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 2 -order 2
402c4762a1bSJed Brown     test:
403c4762a1bSJed Brown       suffix: 4
40430602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,1 -order 2,2
405c4762a1bSJed Brown     test:
406c4762a1bSJed Brown       suffix: 5
40730602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,2 -order 2,3
40830602db0SMatthew G. Knepley 
409c4762a1bSJed Brown   # Spectral ordering 3D 6-11
41030602db0SMatthew G. Knepley   testset:
41130602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,2
41230602db0SMatthew G. Knepley 
413c4762a1bSJed Brown     test:
414c4762a1bSJed Brown       suffix: 6
41530602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 2
416c4762a1bSJed Brown     test:
417c4762a1bSJed Brown       suffix: 7
41830602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 3
419c4762a1bSJed Brown     test:
420c4762a1bSJed Brown       suffix: 8
42130602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 1 -order 5
422c4762a1bSJed Brown     test:
423c4762a1bSJed Brown       suffix: 9
42430602db0SMatthew G. Knepley       args: -num_fields 1 -num_components 2 -order 2
425c4762a1bSJed Brown     test:
426c4762a1bSJed Brown       suffix: 10
42730602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,1 -order 2,2
428c4762a1bSJed Brown     test:
429c4762a1bSJed Brown       suffix: 11
43030602db0SMatthew G. Knepley       args: -num_fields 2 -num_components 1,2 -order 2,3
431c4762a1bSJed Brown 
432c4762a1bSJed Brown TEST*/
433