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 119371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 12c4762a1bSJed Brown PetscInt len; 13c4762a1bSJed Brown PetscBool flg; 14c4762a1bSJed Brown 15c4762a1bSJed Brown PetscFunctionBeginUser; 16c4762a1bSJed Brown options->Nf = 0; 17c4762a1bSJed Brown options->Nc = NULL; 18c4762a1bSJed Brown options->k = NULL; 19c4762a1bSJed Brown 20d0609cedSBarry Smith PetscOptionsBegin(comm, "", "SEM Problem Options", "DMPLEX"); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_fields", "The number of fields", "ex6.c", options->Nf, &options->Nf, NULL, 0)); 22c4762a1bSJed Brown if (options->Nf) { 23c4762a1bSJed Brown len = options->Nf; 249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &options->Nc)); 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-num_components", "The number of components per field", "ex6.c", options->Nc, &len, &flg)); 2663a3b9bcSJacob 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); 27c4762a1bSJed Brown len = options->Nf; 289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &options->k)); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-order", "The spectral order per field", "ex6.c", options->k, &len, &flg)); 3063a3b9bcSJacob 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); 31c4762a1bSJed Brown } 32d0609cedSBarry Smith PetscOptionsEnd(); 33c4762a1bSJed Brown PetscFunctionReturn(0); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown 369371c9d4SSatish Balay static PetscErrorCode LoadData2D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt clSize, Vec u, AppCtx *user) { 37c4762a1bSJed Brown PetscInt i, j, f, c; 38c4762a1bSJed Brown PetscScalar *closure; 39c4762a1bSJed Brown 40c4762a1bSJed Brown PetscFunctionBeginUser; 419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(clSize, &closure)); 42c4762a1bSJed Brown for (j = 0; j < Nj; ++j) { 43c4762a1bSJed Brown for (i = 0; i < Ni; ++i) { 44c4762a1bSJed Brown PetscInt ki, kj, o = 0; 459566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(closure, clSize)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown for (f = 0; f < user->Nf; ++f) { 48c4762a1bSJed Brown PetscInt ioff = i * user->k[f], joff = j * user->k[f]; 49c4762a1bSJed Brown 50c4762a1bSJed Brown for (kj = 0; kj <= user->k[f]; ++kj) { 51c4762a1bSJed Brown for (ki = 0; ki <= user->k[f]; ++ki) { 52*ad540459SPierre Jolivet for (c = 0; c < user->Nc[f]; ++c) closure[o++] = ((kj + joff) * (Ni * user->k[f] + 1) + ki + ioff) * user->Nc[f] + c; 53c4762a1bSJed Brown } 54c4762a1bSJed Brown } 55c4762a1bSJed Brown } 569566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dm, NULL, u, j * Ni + i, closure, INSERT_VALUES)); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown } 599566063dSJacob Faibussowitsch PetscCall(PetscFree(closure)); 60c4762a1bSJed Brown PetscFunctionReturn(0); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 639371c9d4SSatish Balay static PetscErrorCode LoadData3D(DM dm, PetscInt Ni, PetscInt Nj, PetscInt Nk, PetscInt clSize, Vec u, AppCtx *user) { 64c4762a1bSJed Brown PetscInt i, j, k, f, c; 65c4762a1bSJed Brown PetscScalar *closure; 66c4762a1bSJed Brown 67c4762a1bSJed Brown PetscFunctionBeginUser; 689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(clSize, &closure)); 69c4762a1bSJed Brown for (k = 0; k < Nk; ++k) { 70c4762a1bSJed Brown for (j = 0; j < Nj; ++j) { 71c4762a1bSJed Brown for (i = 0; i < Ni; ++i) { 72c4762a1bSJed Brown PetscInt ki, kj, kk, o = 0; 739566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(closure, clSize)); 74c4762a1bSJed Brown 75c4762a1bSJed Brown for (f = 0; f < user->Nf; ++f) { 76c4762a1bSJed Brown PetscInt ioff = i * user->k[f], joff = j * user->k[f], koff = k * user->k[f]; 77c4762a1bSJed Brown 78c4762a1bSJed Brown for (kk = 0; kk <= user->k[f]; ++kk) { 79c4762a1bSJed Brown for (kj = 0; kj <= user->k[f]; ++kj) { 80c4762a1bSJed Brown for (ki = 0; ki <= user->k[f]; ++ki) { 81*ad540459SPierre 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; 82c4762a1bSJed Brown } 83c4762a1bSJed Brown } 84c4762a1bSJed Brown } 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dm, NULL, u, (k * Nj + j) * Ni + i, closure, INSERT_VALUES)); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown } 89c4762a1bSJed Brown } 909566063dSJacob Faibussowitsch PetscCall(PetscFree(closure)); 91c4762a1bSJed Brown PetscFunctionReturn(0); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 949371c9d4SSatish Balay static PetscErrorCode CheckPoint(DM dm, Vec u, PetscInt point, AppCtx *user) { 95c4762a1bSJed Brown PetscSection s; 96c4762a1bSJed Brown PetscScalar *a; 97c4762a1bSJed Brown const PetscScalar *array; 98c4762a1bSJed Brown PetscInt dof, d; 99c4762a1bSJed Brown 100c4762a1bSJed Brown PetscFunctionBeginUser; 1019566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 1029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &array)); 1039566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dm, point, array, &a)); 1049566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(s, point, &dof)); 10563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT ": ", point)); 106c4762a1bSJed Brown for (d = 0; d < dof; ++d) { 1079566063dSJacob Faibussowitsch if (d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 1089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double)PetscRealPart(a[d]))); 109c4762a1bSJed Brown } 1109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 1119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &array)); 112c4762a1bSJed Brown PetscFunctionReturn(0); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 1159371c9d4SSatish Balay static PetscErrorCode ReadData2D(DM dm, Vec u, AppCtx *user) { 116c4762a1bSJed Brown PetscInt cStart, cEnd, cell; 117c4762a1bSJed Brown 118c4762a1bSJed Brown PetscFunctionBeginUser; 1199566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 120c4762a1bSJed Brown for (cell = cStart; cell < cEnd; ++cell) { 121c4762a1bSJed Brown PetscScalar *closure = NULL; 122c4762a1bSJed Brown PetscInt closureSize, ki, kj, f, c, foff = 0; 123c4762a1bSJed Brown 1249566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure)); 12563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT "\n", cell)); 126c4762a1bSJed Brown for (f = 0; f < user->Nf; ++f) { 12763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Field %" PetscInt_FMT "\n", f)); 128c4762a1bSJed Brown for (kj = user->k[f]; kj >= 0; --kj) { 129c4762a1bSJed Brown for (ki = 0; ki <= user->k[f]; ++ki) { 1309566063dSJacob Faibussowitsch if (ki > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " ")); 131c4762a1bSJed Brown for (c = 0; c < user->Nc[f]; ++c) { 1329566063dSJacob Faibussowitsch if (c > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ",")); 1339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%2.0f", (double)PetscRealPart(closure[(kj * (user->k[f] + 1) + ki) * user->Nc[f] + c + foff]))); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown } 1369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 137c4762a1bSJed Brown } 1389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n")); 139c4762a1bSJed Brown foff += PetscSqr(user->k[f] + 1); 140c4762a1bSJed Brown } 1419566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure)); 1429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n")); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown PetscFunctionReturn(0); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 1479371c9d4SSatish Balay static PetscErrorCode ReadData3D(DM dm, Vec u, AppCtx *user) { 148c4762a1bSJed Brown PetscInt cStart, cEnd, cell; 149c4762a1bSJed Brown 150c4762a1bSJed Brown PetscFunctionBeginUser; 1519566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 152c4762a1bSJed Brown for (cell = cStart; cell < cEnd; ++cell) { 153c4762a1bSJed Brown PetscScalar *closure = NULL; 154c4762a1bSJed Brown PetscInt closureSize, ki, kj, kk, f, c, foff = 0; 155c4762a1bSJed Brown 1569566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, u, cell, &closureSize, &closure)); 15763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT "\n", cell)); 158c4762a1bSJed Brown for (f = 0; f < user->Nf; ++f) { 15963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Field %" PetscInt_FMT "\n", f)); 160c4762a1bSJed Brown for (kk = user->k[f]; kk >= 0; --kk) { 161c4762a1bSJed Brown for (kj = user->k[f]; kj >= 0; --kj) { 162c4762a1bSJed Brown for (ki = 0; ki <= user->k[f]; ++ki) { 1639566063dSJacob Faibussowitsch if (ki > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " ")); 164c4762a1bSJed Brown for (c = 0; c < user->Nc[f]; ++c) { 1659566063dSJacob Faibussowitsch if (c > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ",")); 1669566063dSJacob 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]))); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown } 1699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 170c4762a1bSJed Brown } 1719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 172c4762a1bSJed Brown } 1739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n")); 174c4762a1bSJed Brown foff += PetscSqr(user->k[f] + 1); 175c4762a1bSJed Brown } 1769566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, u, cell, &closureSize, &closure)); 1779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n\n")); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown PetscFunctionReturn(0); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 1829371c9d4SSatish Balay static PetscErrorCode SetSymmetries(DM dm, PetscSection s, AppCtx *user) { 18330602db0SMatthew G. Knepley PetscInt dim, f, o, i, j, k, c, d; 184c4762a1bSJed Brown DMLabel depthLabel; 185c4762a1bSJed Brown 186c4762a1bSJed Brown PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1889566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "depth", &depthLabel)); 189c4762a1bSJed Brown for (f = 0; f < user->Nf; f++) { 190c4762a1bSJed Brown PetscSectionSym sym; 191c4762a1bSJed Brown 192c4762a1bSJed 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 */ 1939566063dSJacob Faibussowitsch PetscCall(PetscSectionSymCreateLabel(PetscObjectComm((PetscObject)s), depthLabel, &sym)); 194c4762a1bSJed Brown 19530602db0SMatthew G. Knepley for (d = 0; d <= dim; d++) { 196c4762a1bSJed Brown if (d == 1) { 197c4762a1bSJed Brown PetscInt numDof = user->k[f] - 1; 198c4762a1bSJed Brown PetscInt numComp = user->Nc[f]; 199b5a892a1SMatthew G. Knepley PetscInt minOrnt = -1; 200b5a892a1SMatthew G. Knepley PetscInt maxOrnt = 1; 201c4762a1bSJed Brown PetscInt **perms; 202c4762a1bSJed Brown 2039566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(maxOrnt - minOrnt, &perms)); 204c4762a1bSJed Brown for (o = minOrnt; o < maxOrnt; o++) { 205c4762a1bSJed Brown PetscInt *perm; 206c4762a1bSJed Brown 207b5a892a1SMatthew G. Knepley if (!o) { /* identity */ 208c4762a1bSJed Brown perms[o - minOrnt] = NULL; 209c4762a1bSJed Brown } else { 2109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numDof * numComp, &perm)); 211c4762a1bSJed Brown for (i = numDof - 1, k = 0; i >= 0; i--) { 212c4762a1bSJed Brown for (j = 0; j < numComp; j++, k++) perm[k] = i * numComp + j; 213c4762a1bSJed Brown } 214c4762a1bSJed Brown perms[o - minOrnt] = perm; 215c4762a1bSJed Brown } 216c4762a1bSJed Brown } 2179566063dSJacob Faibussowitsch PetscCall(PetscSectionSymLabelSetStratum(sym, d, numDof * numComp, minOrnt, maxOrnt, PETSC_OWN_POINTER, (const PetscInt **)perms, NULL)); 218c4762a1bSJed Brown } else if (d == 2) { 219c4762a1bSJed Brown PetscInt perEdge = user->k[f] - 1; 220c4762a1bSJed Brown PetscInt numDof = perEdge * perEdge; 221c4762a1bSJed Brown PetscInt numComp = user->Nc[f]; 222c4762a1bSJed Brown PetscInt minOrnt = -4; 223c4762a1bSJed Brown PetscInt maxOrnt = 4; 224c4762a1bSJed Brown PetscInt **perms; 225c4762a1bSJed Brown 2269566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(maxOrnt - minOrnt, &perms)); 227c4762a1bSJed Brown for (o = minOrnt; o < maxOrnt; o++) { 228c4762a1bSJed Brown PetscInt *perm; 229c4762a1bSJed Brown 230c4762a1bSJed Brown if (!o) continue; /* identity */ 2319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numDof * numComp, &perm)); 232c4762a1bSJed Brown /* We want to perm[k] to list which *localArray* position the *sectionArray* position k should go to for the given orientation*/ 233c4762a1bSJed Brown switch (o) { 2349371c9d4SSatish Balay case 0: break; /* identity */ 235b5a892a1SMatthew 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 */ 236c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 237c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 238*ad540459SPierre Jolivet for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * j + i) * numComp + c; 239c4762a1bSJed Brown } 240c4762a1bSJed Brown } 241c4762a1bSJed Brown break; 242b5a892a1SMatthew G. Knepley case -1: /* flip along (-1, 0)--( 1, 0), which swaps edges 0 and 2. This reverses the i variable */ 243c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 244c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 245*ad540459SPierre Jolivet for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + j) * numComp + c; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown } 248c4762a1bSJed Brown break; 249b5a892a1SMatthew 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 */ 250c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 251c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 252*ad540459SPierre Jolivet for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + (perEdge - 1 - i)) * numComp + c; 253c4762a1bSJed Brown } 254c4762a1bSJed Brown } 255c4762a1bSJed Brown break; 256b5a892a1SMatthew G. Knepley case -3: /* flip along ( 0,-1)--( 0, 1), which swaps edges 3 and 1. This reverses the j variable */ 257c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 258c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 259*ad540459SPierre Jolivet for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * i + (perEdge - 1 - j)) * numComp + c; 260c4762a1bSJed Brown } 261c4762a1bSJed Brown } 262c4762a1bSJed Brown break; 263c4762a1bSJed Brown case 1: /* rotate section edge 1 to local edge 0. This swaps the i and j variables and then reverses the j variable */ 264c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 265c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 266*ad540459SPierre Jolivet for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * (perEdge - 1 - j) + i) * numComp + c; 267c4762a1bSJed Brown } 268c4762a1bSJed Brown } 269c4762a1bSJed Brown break; 270c4762a1bSJed Brown case 2: /* rotate section edge 2 to local edge 0. This reverse both i and j variables */ 271c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 272c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 273*ad540459SPierre Jolivet for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * (perEdge - 1 - i) + (perEdge - 1 - j)) * numComp + c; 274c4762a1bSJed Brown } 275c4762a1bSJed Brown } 276c4762a1bSJed Brown break; 277c4762a1bSJed Brown case 3: /* rotate section edge 3 to local edge 0. This swaps the i and j variables and then reverses the i variable */ 278c4762a1bSJed Brown for (i = 0, k = 0; i < perEdge; i++) { 279c4762a1bSJed Brown for (j = 0; j < perEdge; j++, k++) { 280*ad540459SPierre Jolivet for (c = 0; c < numComp; c++) perm[k * numComp + c] = (perEdge * j + (perEdge - 1 - i)) * numComp + c; 281c4762a1bSJed Brown } 282c4762a1bSJed Brown } 283c4762a1bSJed Brown break; 2849371c9d4SSatish Balay default: break; 285c4762a1bSJed Brown } 286c4762a1bSJed Brown perms[o - minOrnt] = perm; 287c4762a1bSJed Brown } 2889566063dSJacob Faibussowitsch PetscCall(PetscSectionSymLabelSetStratum(sym, d, numDof * numComp, minOrnt, maxOrnt, PETSC_OWN_POINTER, (const PetscInt **)perms, NULL)); 289c4762a1bSJed Brown } 290c4762a1bSJed Brown } 2919566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldSym(s, f, sym)); 2929566063dSJacob Faibussowitsch PetscCall(PetscSectionSymDestroy(&sym)); 293c4762a1bSJed Brown } 2949566063dSJacob Faibussowitsch PetscCall(PetscSectionViewFromOptions(s, NULL, "-section_with_sym_view")); 295c4762a1bSJed Brown PetscFunctionReturn(0); 296c4762a1bSJed Brown } 297c4762a1bSJed Brown 2989371c9d4SSatish Balay int main(int argc, char **argv) { 299c4762a1bSJed Brown DM dm; 300c4762a1bSJed Brown PetscSection s; 301c4762a1bSJed Brown Vec u; 302c4762a1bSJed Brown AppCtx user; 30330602db0SMatthew G. Knepley PetscInt dim, size = 0, f; 304c4762a1bSJed Brown 305327415f7SBarry Smith PetscFunctionBeginUser; 3069566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3079566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 3089566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 3099566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 3109566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 3119566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 3129566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 313c4762a1bSJed Brown /* Create a section for SEM order k */ 314c4762a1bSJed Brown { 315c4762a1bSJed Brown PetscInt *numDof, d; 316c4762a1bSJed Brown 3179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user.Nf * (dim + 1), &numDof)); 318c4762a1bSJed Brown for (f = 0; f < user.Nf; ++f) { 31930602db0SMatthew G. Knepley for (d = 0; d <= dim; ++d) numDof[f * (dim + 1) + d] = PetscPowInt(user.k[f] - 1, d) * user.Nc[f]; 320c4762a1bSJed Brown size += PetscPowInt(user.k[f] + 1, d) * user.Nc[f]; 321c4762a1bSJed Brown } 3229566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(dm, user.Nf)); 3239566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(dm, NULL, user.Nc, numDof, 0, NULL, NULL, NULL, NULL, &s)); 3249566063dSJacob Faibussowitsch PetscCall(SetSymmetries(dm, s, &user)); 3259566063dSJacob Faibussowitsch PetscCall(PetscFree(numDof)); 326c4762a1bSJed Brown } 3279566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(dm, s)); 328c4762a1bSJed Brown /* Create spectral ordering and load in data */ 3299566063dSJacob Faibussowitsch PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL)); 3309566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &u)); 33130602db0SMatthew G. Knepley switch (dim) { 3329566063dSJacob Faibussowitsch case 2: PetscCall(LoadData2D(dm, 2, 2, size, u, &user)); break; 3339566063dSJacob Faibussowitsch case 3: PetscCall(LoadData3D(dm, 2, 2, 2, size, u, &user)); break; 334c4762a1bSJed Brown } 335c4762a1bSJed Brown /* Remove ordering and check some values */ 3369566063dSJacob Faibussowitsch PetscCall(PetscSectionSetClosurePermutation(s, (PetscObject)dm, dim, NULL)); 33730602db0SMatthew G. Knepley switch (dim) { 338c4762a1bSJed Brown case 2: 3399566063dSJacob Faibussowitsch PetscCall(CheckPoint(dm, u, 0, &user)); 3409566063dSJacob Faibussowitsch PetscCall(CheckPoint(dm, u, 13, &user)); 3419566063dSJacob Faibussowitsch PetscCall(CheckPoint(dm, u, 15, &user)); 3429566063dSJacob Faibussowitsch PetscCall(CheckPoint(dm, u, 19, &user)); 343c4762a1bSJed Brown break; 344c4762a1bSJed Brown case 3: 3459566063dSJacob Faibussowitsch PetscCall(CheckPoint(dm, u, 0, &user)); 3469566063dSJacob Faibussowitsch PetscCall(CheckPoint(dm, u, 13, &user)); 3479566063dSJacob Faibussowitsch PetscCall(CheckPoint(dm, u, 15, &user)); 3489566063dSJacob Faibussowitsch PetscCall(CheckPoint(dm, u, 19, &user)); 349c4762a1bSJed Brown break; 350c4762a1bSJed Brown } 351c4762a1bSJed Brown /* Recreate spectral ordering and read out data */ 3529566063dSJacob Faibussowitsch PetscCall(DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, s)); 35330602db0SMatthew G. Knepley switch (dim) { 3549566063dSJacob Faibussowitsch case 2: PetscCall(ReadData2D(dm, u, &user)); break; 3559566063dSJacob Faibussowitsch case 3: PetscCall(ReadData3D(dm, u, &user)); break; 356c4762a1bSJed Brown } 3579566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &u)); 3589566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s)); 3599566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3609566063dSJacob Faibussowitsch PetscCall(PetscFree(user.Nc)); 3619566063dSJacob Faibussowitsch PetscCall(PetscFree(user.k)); 3629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 363b122ec5aSJacob Faibussowitsch return 0; 364c4762a1bSJed Brown } 365c4762a1bSJed Brown 366c4762a1bSJed Brown /*TEST 367c4762a1bSJed Brown 368c4762a1bSJed Brown # Spectral ordering 2D 0-5 36930602db0SMatthew G. Knepley testset: 37030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 37130602db0SMatthew G. Knepley 372c4762a1bSJed Brown test: 373c4762a1bSJed Brown suffix: 0 37430602db0SMatthew G. Knepley args: -num_fields 1 -num_components 1 -order 2 375c4762a1bSJed Brown test: 376c4762a1bSJed Brown suffix: 1 37730602db0SMatthew G. Knepley args: -num_fields 1 -num_components 1 -order 3 378c4762a1bSJed Brown test: 379c4762a1bSJed Brown suffix: 2 38030602db0SMatthew G. Knepley args: -num_fields 1 -num_components 1 -order 5 381c4762a1bSJed Brown test: 382c4762a1bSJed Brown suffix: 3 38330602db0SMatthew G. Knepley args: -num_fields 1 -num_components 2 -order 2 384c4762a1bSJed Brown test: 385c4762a1bSJed Brown suffix: 4 38630602db0SMatthew G. Knepley args: -num_fields 2 -num_components 1,1 -order 2,2 387c4762a1bSJed Brown test: 388c4762a1bSJed Brown suffix: 5 38930602db0SMatthew G. Knepley args: -num_fields 2 -num_components 1,2 -order 2,3 39030602db0SMatthew G. Knepley 391c4762a1bSJed Brown # Spectral ordering 3D 6-11 39230602db0SMatthew G. Knepley testset: 39330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 2,2,2 39430602db0SMatthew G. Knepley 395c4762a1bSJed Brown test: 396c4762a1bSJed Brown suffix: 6 39730602db0SMatthew G. Knepley args: -num_fields 1 -num_components 1 -order 2 398c4762a1bSJed Brown test: 399c4762a1bSJed Brown suffix: 7 40030602db0SMatthew G. Knepley args: -num_fields 1 -num_components 1 -order 3 401c4762a1bSJed Brown test: 402c4762a1bSJed Brown suffix: 8 40330602db0SMatthew G. Knepley args: -num_fields 1 -num_components 1 -order 5 404c4762a1bSJed Brown test: 405c4762a1bSJed Brown suffix: 9 40630602db0SMatthew G. Knepley args: -num_fields 1 -num_components 2 -order 2 407c4762a1bSJed Brown test: 408c4762a1bSJed Brown suffix: 10 40930602db0SMatthew G. Knepley args: -num_fields 2 -num_components 1,1 -order 2,2 410c4762a1bSJed Brown test: 411c4762a1bSJed Brown suffix: 11 41230602db0SMatthew G. Knepley args: -num_fields 2 -num_components 1,2 -order 2,3 413c4762a1bSJed Brown 414c4762a1bSJed Brown TEST*/ 415