xref: /petsc/src/dm/tests/ex50.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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) {
42*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm));
43*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetType(dm, DMPLEX));
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetFromOptions(dm));
45c4762a1bSJed Brown   } else {
46*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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   }
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(dm));
49c4762a1bSJed Brown   if (!plex) {
50*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
63*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-order",&dof,NULL));
64c4762a1bSJed Brown     dof += 1;
65c4762a1bSJed Brown 
66*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(dof,&nodes));
67*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(dof,&weights));
68*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTGaussLobattoLegendreQuadrature(dof,PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA,nodes,weights));
69*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinatesLocal(dm,&cv));
70*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinateDM(dm,&cdm));
71c4762a1bSJed Brown     if (plex) {
72c4762a1bSJed Brown       PetscInt cEnd,cStart;
73c4762a1bSJed Brown 
74*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd));
75*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetCoordinateSection(dm,&cSec));
76c4762a1bSJed Brown       nel  = cEnd - cStart;
77c4762a1bSJed Brown       nex  = nel;
78c4762a1bSJed Brown       ney  = 1;
79c4762a1bSJed Brown       nez  = 1;
80c4762a1bSJed Brown     } else {
81*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDAVecGetArray(cdm,cv,&_coords));
82*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDAGetElementsCorners(dm,&gx,&gy,&gz));
83*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDAGetElementsSizes(dm,&nex,&ney,&nez));
84c4762a1bSJed Brown       nel  = nex*ney*nez;
85c4762a1bSJed Brown     }
86*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreate(PETSC_COMM_WORLD,&v));
87*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetSizes(v,3*dof*dof*dof*nel,PETSC_DECIDE));
88*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetType(v,VECSTANDARD));
89*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
97*5f80ce2aSJacob Faibussowitsch             CHKERRQ(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];
101*5f80ce2aSJacob Faibussowitsch             CHKERRQ(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];
115*5f80ce2aSJacob Faibussowitsch                 CHKERRQ(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) {
127*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDAVecRestoreArray(cdm,cv,&_coords));
128c4762a1bSJed Brown     }
129*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(v,&c));
130*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSNPrintf(name,sizeof(name),"FiniteElementCollection: L2_T1_3D_P%D",dof-1));
131*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject)v,name));
132*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)v));
133*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&v));
134*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(nodes));
135*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(weights));
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMViewFromOptions(dm,NULL,"-view"));
137c4762a1bSJed Brown   } else { /* map the whole domain to a sphere */
138*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinates(dm,&v));
139*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetLocalSize(v,&nl));
140*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(v,&c));
141c4762a1bSJed Brown     for (i=0;i<nl/3;i++) {
142*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MapPoint(c+3*i,c+3*i));
143c4762a1bSJed Brown     }
144*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(v,&c));
145*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetCoordinates(dm,v));
146*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMViewFromOptions(dm,NULL,"-view"));
147c4762a1bSJed Brown   }
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
149c4762a1bSJed Brown   PetscFunctionReturn(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown int main(int argc, char *argv[])
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   PetscErrorCode ierr;
155c4762a1bSJed Brown   PetscBool      ho = PETSC_TRUE;
156c4762a1bSJed Brown   PetscBool      plex = PETSC_FALSE;
157c4762a1bSJed Brown   PetscInt       cells[3] = {2,2,2};
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-ho",&ho,NULL));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-plex",&plex,NULL));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nex",&cells[0],NULL));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ney",&cells[1],NULL));
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nez",&cells[2],NULL));
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(test_3d(cells,plex,ho));
166c4762a1bSJed Brown   ierr = PetscFinalize();
167c4762a1bSJed Brown   return ierr;
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