1 static char help[] = "Test GLVis high-order support with DMDAs\n\n"; 2 3 #include <petscdm.h> 4 #include <petscdmda.h> 5 #include <petscdmplex.h> 6 #include <petscdt.h> 7 8 static PetscErrorCode MapPoint(PetscScalar xyz[], PetscScalar mxyz[]) { 9 PetscScalar x, y, z, x2, y2, z2; 10 11 x = xyz[0]; 12 y = xyz[1]; 13 z = xyz[2]; 14 x2 = x * x; 15 y2 = y * y; 16 z2 = z * z; 17 mxyz[0] = x * PetscSqrtScalar(1.0 - y2 / 2.0 - z2 / 2.0 + y2 * z2 / 3.0); 18 mxyz[1] = y * PetscSqrtScalar(1.0 - z2 / 2.0 - x2 / 2.0 + z2 * x2 / 3.0); 19 mxyz[2] = z * PetscSqrtScalar(1.0 - x2 / 2.0 - y2 / 2.0 + x2 * y2 / 3.0); 20 return 0; 21 } 22 23 static PetscErrorCode test_3d(PetscInt cells[], PetscBool plex, PetscBool ho) { 24 DM dm; 25 Vec v; 26 PetscScalar *c; 27 PetscInt nl, i; 28 PetscReal u[3] = {1.0, 1.0, 1.0}, l[3] = {-1.0, -1.0, -1.0}; 29 30 PetscFunctionBeginUser; 31 if (ho) { 32 u[0] = 2. * cells[0]; 33 u[1] = 2. * cells[1]; 34 u[2] = 2. * cells[2]; 35 l[0] = 0.; 36 l[1] = 0.; 37 l[2] = 0.; 38 } 39 if (plex) { 40 PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 41 PetscCall(DMSetType(dm, DMPLEX)); 42 PetscCall(DMSetFromOptions(dm)); 43 } else { 44 PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, cells[0] + 1, cells[1] + 1, cells[2] + 1, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, NULL, &dm)); 45 } 46 PetscCall(DMSetUp(dm)); 47 if (!plex) { PetscCall(DMDASetUniformCoordinates(dm, l[0], u[0], l[1], u[1], l[2], u[2])); } 48 if (ho) { /* each element mapped to a sphere */ 49 DM cdm; 50 Vec cv; 51 PetscSection cSec; 52 DMDACoor3d ***_coords; 53 PetscScalar shift[3], *cptr; 54 PetscInt nel, dof = 3, nex, ney, nez, gx = 0, gy = 0, gz = 0; 55 PetscInt i, j, k, pi, pj, pk; 56 PetscReal *nodes, *weights; 57 char name[256]; 58 59 PetscCall(PetscOptionsGetInt(NULL, NULL, "-order", &dof, NULL)); 60 dof += 1; 61 62 PetscCall(PetscMalloc1(dof, &nodes)); 63 PetscCall(PetscMalloc1(dof, &weights)); 64 PetscCall(PetscDTGaussLobattoLegendreQuadrature(dof, PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, nodes, weights)); 65 PetscCall(DMGetCoordinatesLocal(dm, &cv)); 66 PetscCall(DMGetCoordinateDM(dm, &cdm)); 67 if (plex) { 68 PetscInt cEnd, cStart; 69 70 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 71 PetscCall(DMGetCoordinateSection(dm, &cSec)); 72 nel = cEnd - cStart; 73 nex = nel; 74 ney = 1; 75 nez = 1; 76 } else { 77 PetscCall(DMDAVecGetArray(cdm, cv, &_coords)); 78 PetscCall(DMDAGetElementsCorners(dm, &gx, &gy, &gz)); 79 PetscCall(DMDAGetElementsSizes(dm, &nex, &ney, &nez)); 80 nel = nex * ney * nez; 81 } 82 PetscCall(VecCreate(PETSC_COMM_WORLD, &v)); 83 PetscCall(VecSetSizes(v, 3 * dof * dof * dof * nel, PETSC_DECIDE)); 84 PetscCall(VecSetType(v, VECSTANDARD)); 85 PetscCall(VecGetArray(v, &c)); 86 cptr = c; 87 for (k = gz; k < gz + nez; k++) { 88 for (j = gy; j < gy + ney; j++) { 89 for (i = gx; i < gx + nex; i++) { 90 if (plex) { 91 PetscScalar *t = NULL; 92 93 PetscCall(DMPlexVecGetClosure(dm, cSec, cv, i, NULL, &t)); 94 shift[0] = t[0]; 95 shift[1] = t[1]; 96 shift[2] = t[2]; 97 PetscCall(DMPlexVecRestoreClosure(dm, cSec, cv, i, NULL, &t)); 98 } else { 99 shift[0] = _coords[k][j][i].x; 100 shift[1] = _coords[k][j][i].y; 101 shift[2] = _coords[k][j][i].z; 102 } 103 for (pk = 0; pk < dof; pk++) { 104 PetscScalar xyz[3]; 105 106 xyz[2] = nodes[pk]; 107 for (pj = 0; pj < dof; pj++) { 108 xyz[1] = nodes[pj]; 109 for (pi = 0; pi < dof; pi++) { 110 xyz[0] = nodes[pi]; 111 PetscCall(MapPoint(xyz, cptr)); 112 cptr[0] += shift[0]; 113 cptr[1] += shift[1]; 114 cptr[2] += shift[2]; 115 cptr += 3; 116 } 117 } 118 } 119 } 120 } 121 } 122 if (!plex) { PetscCall(DMDAVecRestoreArray(cdm, cv, &_coords)); } 123 PetscCall(VecRestoreArray(v, &c)); 124 PetscCall(PetscSNPrintf(name, sizeof(name), "FiniteElementCollection: L2_T1_3D_P%" PetscInt_FMT, dof - 1)); 125 PetscCall(PetscObjectSetName((PetscObject)v, name)); 126 PetscCall(PetscObjectCompose((PetscObject)dm, "_glvis_mesh_coords", (PetscObject)v)); 127 PetscCall(VecDestroy(&v)); 128 PetscCall(PetscFree(nodes)); 129 PetscCall(PetscFree(weights)); 130 PetscCall(DMViewFromOptions(dm, NULL, "-view")); 131 } else { /* map the whole domain to a sphere */ 132 PetscCall(DMGetCoordinates(dm, &v)); 133 PetscCall(VecGetLocalSize(v, &nl)); 134 PetscCall(VecGetArray(v, &c)); 135 for (i = 0; i < nl / 3; i++) { PetscCall(MapPoint(c + 3 * i, c + 3 * i)); } 136 PetscCall(VecRestoreArray(v, &c)); 137 PetscCall(DMSetCoordinates(dm, v)); 138 PetscCall(DMViewFromOptions(dm, NULL, "-view")); 139 } 140 PetscCall(DMDestroy(&dm)); 141 PetscFunctionReturn(0); 142 } 143 144 int main(int argc, char *argv[]) { 145 PetscBool ho = PETSC_TRUE; 146 PetscBool plex = PETSC_FALSE; 147 PetscInt cells[3] = {2, 2, 2}; 148 149 PetscFunctionBeginUser; 150 PetscCall(PetscInitialize(&argc, &argv, 0, help)); 151 PetscCall(PetscOptionsGetBool(NULL, NULL, "-ho", &ho, NULL)); 152 PetscCall(PetscOptionsGetBool(NULL, NULL, "-plex", &plex, NULL)); 153 PetscCall(PetscOptionsGetInt(NULL, NULL, "-nex", &cells[0], NULL)); 154 PetscCall(PetscOptionsGetInt(NULL, NULL, "-ney", &cells[1], NULL)); 155 PetscCall(PetscOptionsGetInt(NULL, NULL, "-nez", &cells[2], NULL)); 156 PetscCall(test_3d(cells, plex, ho)); 157 PetscCall(PetscFinalize()); 158 return 0; 159 } 160 161 /*TEST 162 163 testset: 164 nsize: 1 165 args: -view glvis: 166 167 test: 168 suffix: dmda_glvis_ho 169 args: -nex 1 -ney 1 -nez 1 170 171 test: 172 suffix: dmda_glvis_lo 173 args: -ho 0 174 175 test: 176 suffix: dmplex_glvis_ho 177 args: -nex 1 -ney 1 -nez 1 178 179 test: 180 suffix: dmplex_glvis_lo 181 args: -ho 0 182 183 TEST*/ 184