Actual source code: grglvis.c
petsc-3.8.2 2017-11-09
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 if present */
162: if (viewer) {
163: DMDAGetFieldNames(da,(const char * const **)&dafieldname);
164: DMDAGetNumVerticesGhosted(daview,&M,&N,&P);
165: PetscMalloc5(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname);
166: for (i=0;i<dof;i++) bss[i] = 1;
167: nf = dof;
169: PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");
170: PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);
171: PetscOptionsEnd();
172: if (bsset) {
173: PetscInt t;
174: for (i=0,t=0;i<nf;i++) t += bss[i];
175: if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
176: } else nf = dof;
178: for (i=0,s=0;i<nf;i++) {
179: PetscStrallocpy(fec,&fec_type[i]);
180: if (bss[i] == 1) {
181: PetscStrallocpy(dafieldname[s],&fieldname[i]);
182: } else {
183: PetscInt b;
184: size_t tlen = 9; /* "Vector-" + end */
185: for (b=0;b<bss[i];b++) {
186: size_t len;
187: PetscStrlen(dafieldname[s+b],&len);
188: tlen += len + 1; /* field + "-" */
189: }
190: PetscMalloc1(tlen,&fieldname[i]);
191: PetscStrcpy(fieldname[i],"Vector-");
192: for (b=0;b<bss[i]-1;b++) {
193: PetscStrcat(fieldname[i],dafieldname[s+b]);
194: PetscStrcat(fieldname[i],"-");
195: }
196: PetscStrcat(fieldname[i],dafieldname[s+b]);
197: }
198: dims[i] = dim;
199: nlocal[i] = M*N*P*bss[i];
200: s += bss[i];
201: }
203: /* the viewer context takes ownership of xlocal (and the the properly ghosted DMDA associated with it) */
204: PetscNew(&ctx);
205: ctx->xlocal = xlocal;
207: PetscViewerGLVisSetFields(viewer,nf,(const char**)fieldname,(const char**)fec_type,nlocal,bss,dims,DMDASampleGLVisFields_Private,ctx,DMDADestroyGLVisViewerCtx_Private);
208: for (i=0;i<nf;i++) {
209: PetscFree(fec_type[i]);
210: PetscFree(fieldname[i]);
211: }
212: PetscFree5(fec_type,nlocal,bss,dims,fieldname);
213: }
214: DMDestroy(&dacoord);
215: DMDestroy(&daview);
216: return(0);
217: }
219: static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
220: {
221: DM da;
222: Vec xcoorl;
223: PetscMPIInt commsize;
224: const PetscScalar *array;
225: PetscContainer glvis_container;
226: PetscInt dim,sdim,i,vid[8],mid,cid;
227: PetscInt sx,sy,sz,ie,je,ke,ien,jen,ken;
228: PetscInt gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
229: PetscBool enabled = PETSC_TRUE, isascii;
230: PetscErrorCode ierr;
235: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
236: if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
237: MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);
238: if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
239: DMGetDimension(dm,&dim);
241: /* get container: determines if a process visualizes is portion of the data or not */
242: PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);
243: if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
244: {
245: PetscViewerGLVisInfo glvis_info;
246: PetscContainerGetPointer(glvis_container,(void**)&glvis_info);
247: enabled = glvis_info->enabled;
248: }
249: PetscObjectQuery((PetscObject)dm,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);
250: /* this can happen if we are calling DMView outside of VecView_GLVis */
251: if (!xcoorl) {DMSetUpGLVisViewer_DMDA((PetscObject)dm,NULL);}
252: PetscObjectQuery((PetscObject)dm,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);
253: if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
254: VecGetDM(xcoorl,&da);
255: if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
256: DMGetCoordinateDim(da,&sdim);
258: PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");
259: PetscViewerASCIIPrintf(viewer,"\ndimension\n");
260: PetscViewerASCIIPrintf(viewer,"%D\n",dim);
262: if (!enabled) {
263: PetscViewerASCIIPrintf(viewer,"\nelements\n");
264: PetscViewerASCIIPrintf(viewer,"%D\n",0);
265: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
266: PetscViewerASCIIPrintf(viewer,"%D\n",0);
267: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
268: PetscViewerASCIIPrintf(viewer,"%D\n",0);
269: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
270: return(0);
271: }
273: DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
274: i = ien;
275: if (dim > 1) i *= jen;
276: if (dim > 2) i *= ken;
277: PetscViewerASCIIPrintf(viewer,"\nelements\n");
278: PetscViewerASCIIPrintf(viewer,"%D\n",i);
279: switch (dim) {
280: case 1:
281: for (ie = 0; ie < ien; ie++) {
282: vid[0] = ie;
283: vid[1] = ie+1;
284: mid = 1; /* material id */
285: cid = 1; /* segment */
286: PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);
287: }
288: break;
289: case 2:
290: for (je = 0; je < jen; je++) {
291: for (ie = 0; ie < ien; ie++) {
292: vid[0] = je*(ien+1) + ie;
293: vid[1] = je*(ien+1) + ie+1;
294: vid[2] = (je+1)*(ien+1) + ie+1;
295: vid[3] = (je+1)*(ien+1) + ie;
296: mid = 1; /* material id */
297: cid = 3; /* quad */
298: PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);
299: }
300: }
301: break;
302: case 3:
303: for (ke = 0; ke < ken; ke++) {
304: for (je = 0; je < jen; je++) {
305: for (ie = 0; ie < ien; ie++) {
306: vid[0] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie;
307: vid[1] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie+1;
308: vid[2] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
309: vid[3] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
310: vid[4] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie;
311: vid[5] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie+1;
312: vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
313: vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
314: mid = 1; /* material id */
315: cid = 5; /* hex */
316: 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]);
317: }
318: }
319: }
320: break;
321: default:
322: SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
323: break;
324: }
325: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
326: PetscViewerASCIIPrintf(viewer,"%D\n",0);
328: DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);
329: DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);
330: DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);
331: kst = gsz != sz ? 1 : 0;
332: jst = gsy != sy ? 1 : 0;
333: ist = gsx != sx ? 1 : 0;
334: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
335: PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);
336: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
337: VecGetArrayRead(xcoorl,&array);
338: for (ke = kst; ke < kst + ken; ke++) {
339: for (je = jst; je < jst + jen; je++) {
340: for (ie = ist; ie < ist + ien; ie++) {
341: PetscInt d;
343: i = ke * gm * gn + je * gm + ie;
344: for (d=0;d<sdim-1;d++) {
345: PetscViewerASCIIPrintf(viewer,"%g ",PetscRealPart(array[sdim*i+d]));
346: }
347: PetscViewerASCIIPrintf(viewer,"%g\n",PetscRealPart(array[sdim*i+d]));
348: }
349: }
350: }
351: VecRestoreArrayRead(xcoorl,&array);
352: return(0);
353: }
355: /* 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 */
356: PETSC_INTERN PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
357: {
359: PetscBool isglvis,isascii;
364: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);
365: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
366: if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
367: if (isglvis) {
368: PetscViewer view;
369: PetscViewerGLVisType type;
371: PetscViewerGLVisGetType_Private(viewer,&type);
372: PetscViewerGLVisGetDMWindow_Private(viewer,&view);
373: if (view) { /* in the socket case, it may happen that the connection failed */
374: if (type == PETSC_VIEWER_GLVIS_SOCKET) {
375: PetscMPIInt size,rank;
376: MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);
377: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
378: PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);
379: }
380: DMDAView_GLVis_ASCII(dm,view);
381: PetscViewerFlush(view);
382: if (type == PETSC_VIEWER_GLVIS_SOCKET) {
383: PetscInt dim;
384: const char* name;
386: PetscObjectGetName((PetscObject)dm,&name);
387: DMGetDimension(dm,&dim);
388: PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);
389: PetscBarrier((PetscObject)dm);
390: }
391: }
392: } else {
393: DMDAView_GLVis_ASCII(dm,viewer);
394: }
395: return(0);
396: }