Actual source code: grglvis.c
petsc-3.8.0 2017-09-26
1: /* Routines to visualize DMDAs and fields through GLVis */
3: #include <petsc/private/dmdaimpl.h>
4: #include <petsc/private/glvisviewerimpl.h>
6: typedef struct {
7: Vec xlocal;
8: } DMDAGLVisViewerCtx;
10: static PetscErrorCode DMDADestroyGLVisViewerCtx_Private(void *vctx)
11: {
12: DMDAGLVisViewerCtx *ctx = (DMDAGLVisViewerCtx*)vctx;
13: PetscErrorCode ierr;
16: VecDestroy(&ctx->xlocal);
17: PetscFree(vctx);
18: return(0);
19: }
21: /* all but the last proc per dimension claim the ghosted node */
22: static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez)
23: {
24: PetscInt M,N,P,sx,sy,sz,ien,jen,ken;
28: /* Appease -Wmaybe-uninitialized */
29: if (nex) *nex = -1;
30: if (ney) *ney = -1;
31: if (nez) *nez = -1;
32: DMDAGetInfo(da,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
33: DMDAGetCorners(da,&sx,&sy,&sz,&ien,&jen,&ken);
34: if (sx+ien == M) ien--;
35: if (sy+jen == N) jen--;
36: if (sz+ken == P) ken--;
37: if (nex) *nex = ien;
38: if (ney) *ney = jen;
39: if (nez) *nez = ken;
40: return(0);
41: }
43: /* inherits number of vertices from DMDAGetNumElementsGhosted */
44: static PetscErrorCode DMDAGetNumVerticesGhosted(DM da, PetscInt *nvx, PetscInt *nvy, PetscInt *nvz)
45: {
46: PetscInt ien,jen,ken,dim;
50: DMGetDimension(da,&dim);
51: DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
52: ien = ien+1;
53: jen = dim > 1 ? jen+1 : 1;
54: ken = dim > 2 ? ken+1 : 1;
55: if (nvx) *nvx = ien;
56: if (nvy) *nvy = jen;
57: if (nvz) *nvz = ken;
58: return(0);
59: }
61: static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx)
62: {
63: DM da;
64: DMDAGLVisViewerCtx *ctx = (DMDAGLVisViewerCtx*)vctx;
65: const PetscScalar *array;
66: PetscScalar **arrayf;
67: PetscInt i,f,ii,ien,jen,ken,ie,je,ke,bs,*bss;
68: PetscInt sx,sy,sz,gsx,gsy,gsz,ist,jst,kst,gm,gn,gp;
69: PetscErrorCode ierr;
72: VecGetDM(ctx->xlocal,&da);
73: if (!da) SETERRQ(PetscObjectComm(oX),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
74: VecGetBlockSize(ctx->xlocal,&bs);
75: DMGlobalToLocalBegin(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);
76: DMGlobalToLocalEnd(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);
77: DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);
78: DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);
79: DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);
80: kst = gsz != sz ? 1 : 0;
81: jst = gsy != sy ? 1 : 0;
82: ist = gsx != sx ? 1 : 0;
83: PetscMalloc2(nf,&arrayf,nf,&bss);
84: VecGetArrayRead(ctx->xlocal,&array);
85: for (f=0;f<nf;f++) {
86: VecGetBlockSize((Vec)oXf[f],&bss[f]);
87: VecGetArray((Vec)oXf[f],&arrayf[f]);
88: }
89: for (ke = kst, ii = 0; ke < kst + ken; ke++) {
90: for (je = jst; je < jst + jen; je++) {
91: for (ie = ist; ie < ist + ien; ie++) {
92: PetscInt cf,b;
93: i = ke * gm * gn + je * gm + ie;
94: for (f=0,cf=0;f<nf;f++)
95: for (b=0;b<bss[f];b++)
96: arrayf[f][bss[f]*ii+b] = array[i*bs+cf++];
97: ii++;
98: }
99: }
100: }
101: for (f=0;f<nf;f++) { VecRestoreArray((Vec)oXf[f],&arrayf[f]); }
102: VecRestoreArrayRead(ctx->xlocal,&array);
103: PetscFree2(arrayf,bss);
104: return(0);
105: }
107: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
108: {
109: DM da = (DM)oda,daview,dacoord;
110: Vec xcoor,xcoorl,xlocal;
111: DMDAGLVisViewerCtx *ctx;
112: const char **dafieldname;
113: char **fec_type,**fieldname,fec[64];
114: const PetscInt *lx,*ly,*lz;
115: PetscInt *nlocal,*bss,*dims;
116: PetscInt dim,M,N,P,m,n,p,dof,s,i,nf;
117: PetscBool bsset;
118: PetscErrorCode ierr;
121: /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
122: DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL);
123: DMDAGetOwnershipRanges(da,&lx,&ly,&lz);
124: DMGetCoordinates(da,&xcoor);
125: PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n");
126: switch (dim) {
127: case 1:
128: DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview);
129: DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord);
130: PetscStrcpy(fec,"FiniteElementCollection: H1_1D_P1");
131: break;
132: case 2:
133: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,dof,1,lx,ly,&daview);
134: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,2,1,lx,ly,&dacoord);
135: PetscStrcpy(fec,"FiniteElementCollection: H1_2D_P1");
136: break;
137: case 3:
138: DMDACreate3d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,P,m,n,p,dof,1,lx,ly,lz,&daview);
139: DMDACreate3d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,P,m,n,p,3,1,lx,ly,lz,&dacoord);
140: PetscStrcpy(fec,"FiniteElementCollection: H1_3D_P1");
141: break;
142: default:
143: SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
144: break;
145: }
146: DMSetUp(daview);
147: DMSetUp(dacoord);
148: if (!xcoor) {
149: DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);
150: DMGetCoordinates(daview,&xcoor);
151: }
152: DMCreateLocalVector(daview,&xlocal);
153: DMCreateLocalVector(dacoord,&xcoorl);
154: DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);
155: DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);
157: /* xcoorl is composed with the original DMDA, ghosted coordinate DMDA is only available through this vector */
158: PetscObjectCompose(oda,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);
159: PetscObjectDereference((PetscObject)xcoorl);
161: /* customize the viewer */
162: DMDAGetFieldNames(da,(const char * const **)&dafieldname);
163: DMDAGetNumVerticesGhosted(daview,&M,&N,&P);
164: PetscMalloc5(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname);
165: for (i=0;i<dof;i++) bss[i] = 1;
166: nf = dof;
168: PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");
169: PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);
170: PetscOptionsEnd();
171: if (bsset) {
172: PetscInt t;
173: for (i=0,t=0;i<nf;i++) t += bss[i];
174: if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
175: } else nf = dof;
177: for (i=0,s=0;i<nf;i++) {
178: PetscStrallocpy(fec,&fec_type[i]);
179: if (bss[i] == 1) {
180: PetscStrallocpy(dafieldname[s],&fieldname[i]);
181: } else {
182: PetscInt b;
183: size_t tlen = 9; /* "Vector-" + end */
184: for (b=0;b<bss[i];b++) {
185: size_t len;
186: PetscStrlen(dafieldname[s+b],&len);
187: tlen += len + 1; /* field + "-" */
188: }
189: PetscMalloc1(tlen,&fieldname[i]);
190: PetscStrcpy(fieldname[i],"Vector-");
191: for (b=0;b<bss[i]-1;b++) {
192: PetscStrcat(fieldname[i],dafieldname[s+b]);
193: PetscStrcat(fieldname[i],"-");
194: }
195: PetscStrcat(fieldname[i],dafieldname[s+b]);
196: }
197: dims[i] = dim;
198: nlocal[i] = M*N*P*bss[i];
199: s += bss[i];
200: }
202: /* the viewer context takes ownership of xlocal (and the the properly ghosted DMDA associated with it) */
203: PetscNew(&ctx);
204: ctx->xlocal = xlocal;
206: PetscViewerGLVisSetFields(viewer,nf,(const char**)fieldname,(const char**)fec_type,nlocal,bss,dims,DMDASampleGLVisFields_Private,ctx,DMDADestroyGLVisViewerCtx_Private);
207: for (i=0;i<nf;i++) {
208: PetscFree(fec_type[i]);
209: PetscFree(fieldname[i]);
210: }
211: PetscFree5(fec_type,nlocal,bss,dims,fieldname);
212: DMDestroy(&dacoord);
213: DMDestroy(&daview);
214: return(0);
215: }
217: static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
218: {
219: DM da;
220: Vec xcoorl;
221: PetscMPIInt commsize;
222: const PetscScalar *array;
223: PetscContainer glvis_container;
224: PetscInt dim,sdim,i,vid[8],mid,cid;
225: PetscInt sx,sy,sz,ie,je,ke,ien,jen,ken;
226: PetscInt gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
227: PetscBool enabled = PETSC_TRUE, isascii;
228: PetscErrorCode ierr;
233: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
234: if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
235: MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);
236: if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
237: DMGetDimension(dm,&dim);
239: /* get container: determines if a process visualizes is portion of the data or not */
240: PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);
241: if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
242: {
243: PetscViewerGLVisInfo glvis_info;
244: PetscContainerGetPointer(glvis_container,(void**)&glvis_info);
245: enabled = glvis_info->enabled;
246: }
247: PetscObjectQuery((PetscObject)dm,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);
248: if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
249: VecGetDM(xcoorl,&da);
250: if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
251: DMGetCoordinateDim(da,&sdim);
253: PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");
254: PetscViewerASCIIPrintf(viewer,"\ndimension\n");
255: PetscViewerASCIIPrintf(viewer,"%D\n",dim);
257: if (!enabled) {
258: PetscViewerASCIIPrintf(viewer,"\nelements\n");
259: PetscViewerASCIIPrintf(viewer,"%D\n",0);
260: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
261: PetscViewerASCIIPrintf(viewer,"%D\n",0);
262: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
263: PetscViewerASCIIPrintf(viewer,"%D\n",0);
264: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
265: return(0);
266: }
268: DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
269: i = ien;
270: if (dim > 1) i *= jen;
271: if (dim > 2) i *= ken;
272: PetscViewerASCIIPrintf(viewer,"\nelements\n");
273: PetscViewerASCIIPrintf(viewer,"%D\n",i);
274: switch (dim) {
275: case 1:
276: for (ie = 0; ie < ien; ie++) {
277: vid[0] = ie;
278: vid[1] = ie+1;
279: mid = 1; /* material id */
280: cid = 1; /* segment */
281: PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);
282: }
283: break;
284: case 2:
285: for (je = 0; je < jen; je++) {
286: for (ie = 0; ie < ien; ie++) {
287: vid[0] = je*(ien+1) + ie;
288: vid[1] = je*(ien+1) + ie+1;
289: vid[2] = (je+1)*(ien+1) + ie+1;
290: vid[3] = (je+1)*(ien+1) + ie;
291: mid = 1; /* material id */
292: cid = 3; /* quad */
293: PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);
294: }
295: }
296: break;
297: case 3:
298: for (ke = 0; ke < ken; ke++) {
299: for (je = 0; je < jen; je++) {
300: for (ie = 0; ie < ien; ie++) {
301: vid[0] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie;
302: vid[1] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie+1;
303: vid[2] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
304: vid[3] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
305: vid[4] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie;
306: vid[5] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie+1;
307: vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
308: vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
309: mid = 1; /* material id */
310: cid = 5; /* hex */
311: PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3],vid[4],vid[5],vid[6],vid[7]);
312: }
313: }
314: }
315: break;
316: default:
317: SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
318: break;
319: }
320: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
321: PetscViewerASCIIPrintf(viewer,"%D\n",0);
323: DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);
324: DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);
325: DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);
326: kst = gsz != sz ? 1 : 0;
327: jst = gsy != sy ? 1 : 0;
328: ist = gsx != sx ? 1 : 0;
329: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
330: PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);
331: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
332: VecGetArrayRead(xcoorl,&array);
333: for (ke = kst; ke < kst + ken; ke++) {
334: for (je = jst; je < jst + jen; je++) {
335: for (ie = ist; ie < ist + ien; ie++) {
336: PetscInt d;
338: i = ke * gm * gn + je * gm + ie;
339: for (d=0;d<sdim-1;d++) {
340: PetscViewerASCIIPrintf(viewer,"%g ",PetscRealPart(array[sdim*i+d]));
341: }
342: PetscViewerASCIIPrintf(viewer,"%g\n",PetscRealPart(array[sdim*i+d]));
343: }
344: }
345: }
346: VecRestoreArrayRead(xcoorl,&array);
347: return(0);
348: }
350: /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump: duplicated code as in plexglvis.c, should be merged together */
351: PETSC_INTERN PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
352: {
354: PetscBool isglvis,isascii;
359: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);
360: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
361: if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
362: if (isglvis) {
363: PetscViewer view;
364: PetscViewerGLVisType type;
366: PetscViewerGLVisGetType_Private(viewer,&type);
367: PetscViewerGLVisGetDMWindow_Private(viewer,&view);
368: if (view) { /* in the socket case, it may happen that the connection failed */
369: if (type == PETSC_VIEWER_GLVIS_SOCKET) {
370: PetscMPIInt size,rank;
371: MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);
372: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
373: PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);
374: }
375: DMDAView_GLVis_ASCII(dm,view);
376: PetscViewerFlush(view);
377: if (type == PETSC_VIEWER_GLVIS_SOCKET) {
378: PetscInt dim;
379: const char* name;
381: PetscObjectGetName((PetscObject)dm,&name);
382: DMGetDimension(dm,&dim);
383: PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);
384: PetscBarrier((PetscObject)dm);
385: }
386: }
387: } else {
388: DMDAView_GLVis_ASCII(dm,viewer);
389: }
390: return(0);
391: }