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