1c4762a1bSJed Brown static char help[] = "Test GLVis high-order support with DMDAs\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdm.h> 4c4762a1bSJed Brown #include <petscdmda.h> 5c4762a1bSJed Brown #include <petscdmplex.h> 6c4762a1bSJed Brown #include <petscdt.h> 7c4762a1bSJed Brown 89371c9d4SSatish Balay static PetscErrorCode MapPoint(PetscScalar xyz[], PetscScalar mxyz[]) { 9c4762a1bSJed Brown PetscScalar x, y, z, x2, y2, z2; 10c4762a1bSJed Brown 11c4762a1bSJed Brown x = xyz[0]; 12c4762a1bSJed Brown y = xyz[1]; 13c4762a1bSJed Brown z = xyz[2]; 14c4762a1bSJed Brown x2 = x * x; 15c4762a1bSJed Brown y2 = y * y; 16c4762a1bSJed Brown z2 = z * z; 17c4762a1bSJed Brown mxyz[0] = x * PetscSqrtScalar(1.0 - y2 / 2.0 - z2 / 2.0 + y2 * z2 / 3.0); 18c4762a1bSJed Brown mxyz[1] = y * PetscSqrtScalar(1.0 - z2 / 2.0 - x2 / 2.0 + z2 * x2 / 3.0); 19c4762a1bSJed Brown mxyz[2] = z * PetscSqrtScalar(1.0 - x2 / 2.0 - y2 / 2.0 + x2 * y2 / 3.0); 20c4762a1bSJed Brown return 0; 21c4762a1bSJed Brown } 22c4762a1bSJed Brown 239371c9d4SSatish Balay static PetscErrorCode test_3d(PetscInt cells[], PetscBool plex, PetscBool ho) { 24c4762a1bSJed Brown DM dm; 25c4762a1bSJed Brown Vec v; 26c4762a1bSJed Brown PetscScalar *c; 27c4762a1bSJed Brown PetscInt nl, i; 28c4762a1bSJed Brown PetscReal u[3] = {1.0, 1.0, 1.0}, l[3] = {-1.0, -1.0, -1.0}; 29c4762a1bSJed Brown 30c4762a1bSJed Brown PetscFunctionBeginUser; 31c4762a1bSJed Brown if (ho) { 32c4762a1bSJed Brown u[0] = 2. * cells[0]; 33c4762a1bSJed Brown u[1] = 2. * cells[1]; 34c4762a1bSJed Brown u[2] = 2. * cells[2]; 35c4762a1bSJed Brown l[0] = 0.; 36c4762a1bSJed Brown l[1] = 0.; 37c4762a1bSJed Brown l[2] = 0.; 38c4762a1bSJed Brown } 39c4762a1bSJed Brown if (plex) { 409566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 419566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 429566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 43c4762a1bSJed Brown } else { 449566063dSJacob Faibussowitsch 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)); 45c4762a1bSJed Brown } 469566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 47*48a46eb9SPierre Jolivet if (!plex) PetscCall(DMDASetUniformCoordinates(dm, l[0], u[0], l[1], u[1], l[2], u[2])); 48c4762a1bSJed Brown if (ho) { /* each element mapped to a sphere */ 49c4762a1bSJed Brown DM cdm; 50c4762a1bSJed Brown Vec cv; 51c4762a1bSJed Brown PetscSection cSec; 52c4762a1bSJed Brown DMDACoor3d ***_coords; 53c4762a1bSJed Brown PetscScalar shift[3], *cptr; 54c4762a1bSJed Brown PetscInt nel, dof = 3, nex, ney, nez, gx = 0, gy = 0, gz = 0; 55c4762a1bSJed Brown PetscInt i, j, k, pi, pj, pk; 56c4762a1bSJed Brown PetscReal *nodes, *weights; 57c4762a1bSJed Brown char name[256]; 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-order", &dof, NULL)); 60c4762a1bSJed Brown dof += 1; 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, &nodes)); 639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, &weights)); 649566063dSJacob Faibussowitsch PetscCall(PetscDTGaussLobattoLegendreQuadrature(dof, PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA, nodes, weights)); 659566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &cv)); 669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 67c4762a1bSJed Brown if (plex) { 68c4762a1bSJed Brown PetscInt cEnd, cStart; 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 719566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &cSec)); 72c4762a1bSJed Brown nel = cEnd - cStart; 73c4762a1bSJed Brown nex = nel; 74c4762a1bSJed Brown ney = 1; 75c4762a1bSJed Brown nez = 1; 76c4762a1bSJed Brown } else { 779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cdm, cv, &_coords)); 789566063dSJacob Faibussowitsch PetscCall(DMDAGetElementsCorners(dm, &gx, &gy, &gz)); 799566063dSJacob Faibussowitsch PetscCall(DMDAGetElementsSizes(dm, &nex, &ney, &nez)); 80c4762a1bSJed Brown nel = nex * ney * nez; 81c4762a1bSJed Brown } 829566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &v)); 839566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v, 3 * dof * dof * dof * nel, PETSC_DECIDE)); 849566063dSJacob Faibussowitsch PetscCall(VecSetType(v, VECSTANDARD)); 859566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &c)); 86c4762a1bSJed Brown cptr = c; 87c4762a1bSJed Brown for (k = gz; k < gz + nez; k++) { 88c4762a1bSJed Brown for (j = gy; j < gy + ney; j++) { 89c4762a1bSJed Brown for (i = gx; i < gx + nex; i++) { 90c4762a1bSJed Brown if (plex) { 91c4762a1bSJed Brown PetscScalar *t = NULL; 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, cSec, cv, i, NULL, &t)); 94c4762a1bSJed Brown shift[0] = t[0]; 95c4762a1bSJed Brown shift[1] = t[1]; 96c4762a1bSJed Brown shift[2] = t[2]; 979566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, cSec, cv, i, NULL, &t)); 98c4762a1bSJed Brown } else { 99c4762a1bSJed Brown shift[0] = _coords[k][j][i].x; 100c4762a1bSJed Brown shift[1] = _coords[k][j][i].y; 101c4762a1bSJed Brown shift[2] = _coords[k][j][i].z; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown for (pk = 0; pk < dof; pk++) { 104c4762a1bSJed Brown PetscScalar xyz[3]; 105c4762a1bSJed Brown 106c4762a1bSJed Brown xyz[2] = nodes[pk]; 107c4762a1bSJed Brown for (pj = 0; pj < dof; pj++) { 108c4762a1bSJed Brown xyz[1] = nodes[pj]; 109c4762a1bSJed Brown for (pi = 0; pi < dof; pi++) { 110c4762a1bSJed Brown xyz[0] = nodes[pi]; 1119566063dSJacob Faibussowitsch PetscCall(MapPoint(xyz, cptr)); 112c4762a1bSJed Brown cptr[0] += shift[0]; 113c4762a1bSJed Brown cptr[1] += shift[1]; 114c4762a1bSJed Brown cptr[2] += shift[2]; 115c4762a1bSJed Brown cptr += 3; 116c4762a1bSJed Brown } 117c4762a1bSJed Brown } 118c4762a1bSJed Brown } 119c4762a1bSJed Brown } 120c4762a1bSJed Brown } 121c4762a1bSJed Brown } 122*48a46eb9SPierre Jolivet if (!plex) PetscCall(DMDAVecRestoreArray(cdm, cv, &_coords)); 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &c)); 12463a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, sizeof(name), "FiniteElementCollection: L2_T1_3D_P%" PetscInt_FMT, dof - 1)); 1259566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)v, name)); 1269566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "_glvis_mesh_coords", (PetscObject)v)); 1279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1289566063dSJacob Faibussowitsch PetscCall(PetscFree(nodes)); 1299566063dSJacob Faibussowitsch PetscCall(PetscFree(weights)); 1309566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-view")); 131c4762a1bSJed Brown } else { /* map the whole domain to a sphere */ 1329566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &v)); 1339566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &nl)); 1349566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &c)); 135*48a46eb9SPierre Jolivet for (i = 0; i < nl / 3; i++) PetscCall(MapPoint(c + 3 * i, c + 3 * i)); 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &c)); 1379566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(dm, v)); 1389566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-view")); 139c4762a1bSJed Brown } 1409566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 141c4762a1bSJed Brown PetscFunctionReturn(0); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 1449371c9d4SSatish Balay int main(int argc, char *argv[]) { 145c4762a1bSJed Brown PetscBool ho = PETSC_TRUE; 146c4762a1bSJed Brown PetscBool plex = PETSC_FALSE; 147c4762a1bSJed Brown PetscInt cells[3] = {2, 2, 2}; 148c4762a1bSJed Brown 149327415f7SBarry Smith PetscFunctionBeginUser; 1509566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 1519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-ho", &ho, NULL)); 1529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-plex", &plex, NULL)); 1539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nex", &cells[0], NULL)); 1549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ney", &cells[1], NULL)); 1559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nez", &cells[2], NULL)); 1569566063dSJacob Faibussowitsch PetscCall(test_3d(cells, plex, ho)); 1579566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 158b122ec5aSJacob Faibussowitsch return 0; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown /*TEST 162c4762a1bSJed Brown 163c4762a1bSJed Brown testset: 164c4762a1bSJed Brown nsize: 1 165c4762a1bSJed Brown args: -view glvis: 166c4762a1bSJed Brown 167c4762a1bSJed Brown test: 168c4762a1bSJed Brown suffix: dmda_glvis_ho 169c4762a1bSJed Brown args: -nex 1 -ney 1 -nez 1 170c4762a1bSJed Brown 171c4762a1bSJed Brown test: 172c4762a1bSJed Brown suffix: dmda_glvis_lo 173c4762a1bSJed Brown args: -ho 0 174c4762a1bSJed Brown 175c4762a1bSJed Brown test: 176c4762a1bSJed Brown suffix: dmplex_glvis_ho 177c4762a1bSJed Brown args: -nex 1 -ney 1 -nez 1 178c4762a1bSJed Brown 179c4762a1bSJed Brown test: 180c4762a1bSJed Brown suffix: dmplex_glvis_lo 181c4762a1bSJed Brown args: -ho 0 182c4762a1bSJed Brown 183c4762a1bSJed Brown TEST*/ 184