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 8c4762a1bSJed Brown static PetscErrorCode MapPoint(PetscScalar xyz[],PetscScalar mxyz[]) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown PetscScalar x,y,z,x2,y2,z2; 11c4762a1bSJed Brown 12c4762a1bSJed Brown x = xyz[0]; 13c4762a1bSJed Brown y = xyz[1]; 14c4762a1bSJed Brown z = xyz[2]; 15c4762a1bSJed Brown x2 = x*x; 16c4762a1bSJed Brown y2 = y*y; 17c4762a1bSJed Brown z2 = z*z; 18c4762a1bSJed Brown mxyz[0] = x*PetscSqrtScalar(1.0-y2/2.0-z2/2.0 + y2*z2/3.0); 19c4762a1bSJed Brown mxyz[1] = y*PetscSqrtScalar(1.0-z2/2.0-x2/2.0 + z2*x2/3.0); 20c4762a1bSJed Brown mxyz[2] = z*PetscSqrtScalar(1.0-x2/2.0-y2/2.0 + x2*y2/3.0); 21c4762a1bSJed Brown return 0; 22c4762a1bSJed Brown } 23c4762a1bSJed Brown 24c4762a1bSJed Brown static PetscErrorCode test_3d(PetscInt cells[], PetscBool plex, PetscBool ho) 25c4762a1bSJed Brown { 26c4762a1bSJed Brown DM dm; 27c4762a1bSJed Brown Vec v; 28c4762a1bSJed Brown PetscScalar *c; 29c4762a1bSJed Brown PetscInt nl,i; 30c4762a1bSJed Brown PetscReal u[3] = {1.0,1.0,1.0}, l[3] = {-1.0,-1.0,-1.0}; 31c4762a1bSJed Brown 32c4762a1bSJed Brown PetscFunctionBeginUser; 33c4762a1bSJed Brown if (ho) { 34c4762a1bSJed Brown u[0] = 2.*cells[0]; 35c4762a1bSJed Brown u[1] = 2.*cells[1]; 36c4762a1bSJed Brown u[2] = 2.*cells[2]; 37c4762a1bSJed Brown l[0] = 0.; 38c4762a1bSJed Brown l[1] = 0.; 39c4762a1bSJed Brown l[2] = 0.; 40c4762a1bSJed Brown } 41c4762a1bSJed Brown if (plex) { 429566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 439566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 449566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 45c4762a1bSJed Brown } else { 469566063dSJacob 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)); 47c4762a1bSJed Brown } 489566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 49c4762a1bSJed Brown if (!plex) { 509566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dm,l[0],u[0],l[1],u[1],l[2],u[2])); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown if (ho) { /* each element mapped to a sphere */ 53c4762a1bSJed Brown DM cdm; 54c4762a1bSJed Brown Vec cv; 55c4762a1bSJed Brown PetscSection cSec; 56c4762a1bSJed Brown DMDACoor3d ***_coords; 57c4762a1bSJed Brown PetscScalar shift[3],*cptr; 58c4762a1bSJed Brown PetscInt nel,dof = 3,nex,ney,nez,gx = 0,gy = 0,gz = 0; 59c4762a1bSJed Brown PetscInt i,j,k,pi,pj,pk; 60c4762a1bSJed Brown PetscReal *nodes,*weights; 61c4762a1bSJed Brown char name[256]; 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-order",&dof,NULL)); 64c4762a1bSJed Brown dof += 1; 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof,&nodes)); 679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof,&weights)); 689566063dSJacob Faibussowitsch PetscCall(PetscDTGaussLobattoLegendreQuadrature(dof,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,nodes,weights)); 699566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm,&cv)); 709566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm,&cdm)); 71c4762a1bSJed Brown if (plex) { 72c4762a1bSJed Brown PetscInt cEnd,cStart; 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd)); 759566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm,&cSec)); 76c4762a1bSJed Brown nel = cEnd - cStart; 77c4762a1bSJed Brown nex = nel; 78c4762a1bSJed Brown ney = 1; 79c4762a1bSJed Brown nez = 1; 80c4762a1bSJed Brown } else { 819566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cdm,cv,&_coords)); 829566063dSJacob Faibussowitsch PetscCall(DMDAGetElementsCorners(dm,&gx,&gy,&gz)); 839566063dSJacob Faibussowitsch PetscCall(DMDAGetElementsSizes(dm,&nex,&ney,&nez)); 84c4762a1bSJed Brown nel = nex*ney*nez; 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&v)); 879566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v,3*dof*dof*dof*nel,PETSC_DECIDE)); 889566063dSJacob Faibussowitsch PetscCall(VecSetType(v,VECSTANDARD)); 899566063dSJacob Faibussowitsch PetscCall(VecGetArray(v,&c)); 90c4762a1bSJed Brown cptr = c; 91c4762a1bSJed Brown for (k=gz;k<gz+nez;k++) { 92c4762a1bSJed Brown for (j=gy;j<gy+ney;j++) { 93c4762a1bSJed Brown for (i=gx;i<gx+nex;i++) { 94c4762a1bSJed Brown if (plex) { 95c4762a1bSJed Brown PetscScalar *t = NULL; 96c4762a1bSJed Brown 979566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm,cSec,cv,i,NULL,&t)); 98c4762a1bSJed Brown shift[0] = t[0]; 99c4762a1bSJed Brown shift[1] = t[1]; 100c4762a1bSJed Brown shift[2] = t[2]; 1019566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm,cSec,cv,i,NULL,&t)); 102c4762a1bSJed Brown } else { 103c4762a1bSJed Brown shift[0] = _coords[k][j][i].x; 104c4762a1bSJed Brown shift[1] = _coords[k][j][i].y; 105c4762a1bSJed Brown shift[2] = _coords[k][j][i].z; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown for (pk=0;pk<dof;pk++) { 108c4762a1bSJed Brown PetscScalar xyz[3]; 109c4762a1bSJed Brown 110c4762a1bSJed Brown xyz[2] = nodes[pk]; 111c4762a1bSJed Brown for (pj=0;pj<dof;pj++) { 112c4762a1bSJed Brown xyz[1] = nodes[pj]; 113c4762a1bSJed Brown for (pi=0;pi<dof;pi++) { 114c4762a1bSJed Brown xyz[0] = nodes[pi]; 1159566063dSJacob Faibussowitsch PetscCall(MapPoint(xyz,cptr)); 116c4762a1bSJed Brown cptr[0] += shift[0]; 117c4762a1bSJed Brown cptr[1] += shift[1]; 118c4762a1bSJed Brown cptr[2] += shift[2]; 119c4762a1bSJed Brown cptr += 3; 120c4762a1bSJed Brown } 121c4762a1bSJed Brown } 122c4762a1bSJed Brown } 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 125c4762a1bSJed Brown } 126c4762a1bSJed Brown if (!plex) { 1279566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cdm,cv,&_coords)); 128c4762a1bSJed Brown } 1299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v,&c)); 13063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name,sizeof(name),"FiniteElementCollection: L2_T1_3D_P%" PetscInt_FMT,dof-1)); 1319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)v,name)); 1329566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)v)); 1339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1349566063dSJacob Faibussowitsch PetscCall(PetscFree(nodes)); 1359566063dSJacob Faibussowitsch PetscCall(PetscFree(weights)); 1369566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm,NULL,"-view")); 137c4762a1bSJed Brown } else { /* map the whole domain to a sphere */ 1389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm,&v)); 1399566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v,&nl)); 1409566063dSJacob Faibussowitsch PetscCall(VecGetArray(v,&c)); 141c4762a1bSJed Brown for (i=0;i<nl/3;i++) { 1429566063dSJacob Faibussowitsch PetscCall(MapPoint(c+3*i,c+3*i)); 143c4762a1bSJed Brown } 1449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v,&c)); 1459566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(dm,v)); 1469566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm,NULL,"-view")); 147c4762a1bSJed Brown } 1489566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 149c4762a1bSJed Brown PetscFunctionReturn(0); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown int main(int argc, char *argv[]) 153c4762a1bSJed Brown { 154c4762a1bSJed Brown PetscBool ho = PETSC_TRUE; 155c4762a1bSJed Brown PetscBool plex = PETSC_FALSE; 156c4762a1bSJed Brown PetscInt cells[3] = {2,2,2}; 157c4762a1bSJed Brown 158*327415f7SBarry Smith PetscFunctionBeginUser; 1599566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,0,help)); 1609566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-ho",&ho,NULL)); 1619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-plex",&plex,NULL)); 1629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nex",&cells[0],NULL)); 1639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-ney",&cells[1],NULL)); 1649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nez",&cells[2],NULL)); 1659566063dSJacob Faibussowitsch PetscCall(test_3d(cells,plex,ho)); 1669566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 167b122ec5aSJacob Faibussowitsch return 0; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown /*TEST 171c4762a1bSJed Brown 172c4762a1bSJed Brown testset: 173c4762a1bSJed Brown nsize: 1 174c4762a1bSJed Brown args: -view glvis: 175c4762a1bSJed Brown 176c4762a1bSJed Brown test: 177c4762a1bSJed Brown suffix: dmda_glvis_ho 178c4762a1bSJed Brown args: -nex 1 -ney 1 -nez 1 179c4762a1bSJed Brown 180c4762a1bSJed Brown test: 181c4762a1bSJed Brown suffix: dmda_glvis_lo 182c4762a1bSJed Brown args: -ho 0 183c4762a1bSJed Brown 184c4762a1bSJed Brown test: 185c4762a1bSJed Brown suffix: dmplex_glvis_ho 186c4762a1bSJed Brown args: -nex 1 -ney 1 -nez 1 187c4762a1bSJed Brown 188c4762a1bSJed Brown test: 189c4762a1bSJed Brown suffix: dmplex_glvis_lo 190c4762a1bSJed Brown args: -ho 0 191c4762a1bSJed Brown 192c4762a1bSJed Brown TEST*/ 193