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