Actual source code: grglvis.c
petsc-3.9.4 2018-09-11
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: PetscBool ll;
8: } DMDAGhostedGLVisViewerCtx;
10: static PetscErrorCode DMDAGhostedDestroyGLVisViewerCtx_Private(void **vctx)
11: {
15: PetscFree(*vctx);
16: return(0);
17: }
19: typedef struct {
20: Vec xlocal;
21: } DMDAFieldGLVisViewerCtx;
23: static PetscErrorCode DMDAFieldDestroyGLVisViewerCtx_Private(void *vctx)
24: {
25: DMDAFieldGLVisViewerCtx *ctx = (DMDAFieldGLVisViewerCtx*)vctx;
26: PetscErrorCode ierr;
29: VecDestroy(&ctx->xlocal);
30: PetscFree(vctx);
31: return(0);
32: }
34: /*
35: dactx->ll is false -> all but the last proc per dimension claim the ghosted node on the right
36: dactx->ll is true -> all but the first proc per dimension claim the ghosted node on the left
37: */
38: static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez)
39: {
40: DMDAGhostedGLVisViewerCtx *dactx;
41: PetscInt sx,sy,sz,ien,jen,ken;
42: PetscErrorCode ierr;
45: /* Appease -Wmaybe-uninitialized */
46: if (nex) *nex = -1;
47: if (ney) *ney = -1;
48: if (nez) *nez = -1;
49: DMDAGetCorners(da,&sx,&sy,&sz,&ien,&jen,&ken);
50: DMGetApplicationContext(da,(void**)&dactx);
51: if (dactx->ll) {
52: PetscInt dim;
54: DMGetDimension(da,&dim);
55: if (!sx) ien--;
56: if (!sy && dim > 1) jen--;
57: if (!sz && dim > 2) ken--;
58: } else {
59: PetscInt M,N,P;
61: DMDAGetInfo(da,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
62: if (sx+ien == M) ien--;
63: if (sy+jen == N) jen--;
64: if (sz+ken == P) ken--;
65: }
66: if (nex) *nex = ien;
67: if (ney) *ney = jen;
68: if (nez) *nez = ken;
69: return(0);
70: }
72: /* inherits number of vertices from DMDAGetNumElementsGhosted */
73: static PetscErrorCode DMDAGetNumVerticesGhosted(DM da, PetscInt *nvx, PetscInt *nvy, PetscInt *nvz)
74: {
75: PetscInt ien,jen,ken,dim;
79: DMGetDimension(da,&dim);
80: DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
81: ien = ien+1;
82: jen = dim > 1 ? jen+1 : 1;
83: ken = dim > 2 ? ken+1 : 1;
84: if (nvx) *nvx = ien;
85: if (nvy) *nvy = jen;
86: if (nvz) *nvz = ken;
87: return(0);
88: }
90: static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx)
91: {
92: DM da;
93: DMDAFieldGLVisViewerCtx *ctx = (DMDAFieldGLVisViewerCtx*)vctx;
94: DMDAGhostedGLVisViewerCtx *dactx;
95: const PetscScalar *array;
96: PetscScalar **arrayf;
97: PetscInt i,f,ii,ien,jen,ken,ie,je,ke,bs,*bss;
98: PetscInt sx,sy,sz,gsx,gsy,gsz,ist,jst,kst,gm,gn,gp;
99: PetscErrorCode ierr;
102: VecGetDM(ctx->xlocal,&da);
103: if (!da) SETERRQ(PetscObjectComm(oX),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
104: DMGetApplicationContext(da,(void**)&dactx);
105: VecGetBlockSize(ctx->xlocal,&bs);
106: DMGlobalToLocalBegin(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);
107: DMGlobalToLocalEnd(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);
108: DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);
109: DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);
110: DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);
111: if (dactx->ll) {
112: kst = jst = ist = 0;
113: } else {
114: kst = gsz != sz ? 1 : 0;
115: jst = gsy != sy ? 1 : 0;
116: ist = gsx != sx ? 1 : 0;
117: }
118: PetscMalloc2(nf,&arrayf,nf,&bss);
119: VecGetArrayRead(ctx->xlocal,&array);
120: for (f=0;f<nf;f++) {
121: VecGetBlockSize((Vec)oXf[f],&bss[f]);
122: VecGetArray((Vec)oXf[f],&arrayf[f]);
123: }
124: for (ke = kst, ii = 0; ke < kst + ken; ke++) {
125: for (je = jst; je < jst + jen; je++) {
126: for (ie = ist; ie < ist + ien; ie++) {
127: PetscInt cf,b;
128: i = ke * gm * gn + je * gm + ie;
129: for (f=0,cf=0;f<nf;f++)
130: for (b=0;b<bss[f];b++)
131: arrayf[f][bss[f]*ii+b] = array[i*bs+cf++];
132: ii++;
133: }
134: }
135: }
136: for (f=0;f<nf;f++) { VecRestoreArray((Vec)oXf[f],&arrayf[f]); }
137: VecRestoreArrayRead(ctx->xlocal,&array);
138: PetscFree2(arrayf,bss);
139: return(0);
140: }
142: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
143: {
144: DM da = (DM)oda,daview;
148: PetscObjectQuery(oda,"GLVisGraphicsDMDAGhosted",(PetscObject*)&daview);
149: if (!daview) {
150: DMDAGhostedGLVisViewerCtx *dactx;
151: DM dacoord = NULL;
152: Vec xcoor,xcoorl;
153: PetscBool hashocoord = PETSC_FALSE;
154: const PetscInt *lx,*ly,*lz;
155: PetscInt dim,M,N,P,m,n,p,dof,s,i;
157: PetscNew(&dactx);
158: dactx->ll = PETSC_TRUE; /* default to match elements layout obtained by DMDAGetElements */
159: PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");
160: PetscOptionsBool("-viewer_glvis_dm_da_ll","Left-looking subdomain view",NULL,dactx->ll,&dactx->ll,NULL);
161: PetscOptionsEnd();
162: /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
163: DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL);
164: DMDAGetOwnershipRanges(da,&lx,&ly,&lz);
165: PetscObjectQuery((PetscObject)da,"_glvis_mesh_coords",(PetscObject*)&xcoor);
166: if (!xcoor) {
167: DMGetCoordinates(da,&xcoor);
168: } else {
169: hashocoord = PETSC_TRUE;
170: }
171: PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n");
172: switch (dim) {
173: case 1:
174: DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview);
175: if (!hashocoord) {
176: DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord);
177: }
178: break;
179: case 2:
180: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,dof,1,lx,ly,&daview);
181: if (!hashocoord) {
182: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,2,1,lx,ly,&dacoord);
183: }
184: break;
185: case 3:
186: 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);
187: if (!hashocoord) {
188: 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);
189: }
190: break;
191: default:
192: SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
193: break;
194: }
195: DMSetApplicationContext(daview,dactx);
196: DMSetApplicationContextDestroy(daview,DMDAGhostedDestroyGLVisViewerCtx_Private);
197: DMSetUp(daview);
198: if (!xcoor) {
199: DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);
200: DMGetCoordinates(daview,&xcoor);
201: }
202: if (dacoord) {
203: DMSetUp(dacoord);
204: DMCreateLocalVector(dacoord,&xcoorl);
205: DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);
206: DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);
207: DMDestroy(&dacoord);
208: } else {
209: PetscInt ien,jen,ken,nc,nl,cdof,deg;
210: char fecmesh[64];
212: DMDAGetNumElementsGhosted(daview,&ien,&jen,&ken);
213: nc = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
215: VecGetLocalSize(xcoor,&nl);
216: if (nc && nl % nc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible local coordinate size %D and number of cells %D",nl,nc);
217: VecDuplicate(xcoor,&xcoorl);
218: VecCopy(xcoor,xcoorl);
219: VecSetDM(xcoorl,NULL);
220: cdof = nl/(nc*dim);
221: deg = 1;
222: while (1) {
223: PetscInt degd = 1;
224: for (i=0;i<dim;i++) degd *= (deg+1);
225: if (degd > cdof) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell dofs %D",cdof);
226: if (degd == cdof) break;
227: deg++;
228: }
229: PetscSNPrintf(fecmesh,sizeof(fecmesh),"FiniteElementCollection: L2_T1_%DD_P%D",dim,deg);
230: PetscObjectSetName((PetscObject)xcoorl,fecmesh);
231: }
233: /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */
234: PetscObjectCompose((PetscObject)daview,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);
235: PetscObjectDereference((PetscObject)xcoorl);
237: /* daview is composed with the original DMDA */
238: PetscObjectCompose(oda,"GLVisGraphicsDMDAGhosted",(PetscObject)daview);
239: PetscObjectDereference((PetscObject)daview);
240: }
242: /* customize the viewer if present */
243: if (viewer) {
244: DMDAFieldGLVisViewerCtx *ctx;
245: DMDAGhostedGLVisViewerCtx *dactx;
246: char fec[64];
247: Vec xlocal,*Ufield;
248: const char **dafieldname;
249: char **fec_type,**fieldname;
250: PetscInt *nlocal,*bss,*dims;
251: PetscInt dim,M,N,P,dof,s,i,nf;
252: PetscBool bsset;
254: DMGetApplicationContext(daview,(void**)&dactx);
255: DMCreateLocalVector(daview,&xlocal);
256: DMDAGetFieldNames(da,(const char * const **)&dafieldname);
257: DMDAGetNumVerticesGhosted(daview,&M,&N,&P);
258: DMDAGetInfo(daview,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
259: PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: H1_%DD_P1",dim);
260: PetscMalloc6(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname,dof,&Ufield);
261: for (i=0;i<dof;i++) bss[i] = 1;
262: nf = dof;
264: PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Field options","PetscViewer");
265: PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);
266: PetscOptionsEnd();
267: if (bsset) {
268: PetscInt t;
269: for (i=0,t=0;i<nf;i++) t += bss[i];
270: if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
271: } else nf = dof;
273: for (i=0,s=0;i<nf;i++) {
274: PetscStrallocpy(fec,&fec_type[i]);
275: if (bss[i] == 1) {
276: PetscStrallocpy(dafieldname[s],&fieldname[i]);
277: } else {
278: PetscInt b;
279: size_t tlen = 9; /* "Vector-" + end */
280: for (b=0;b<bss[i];b++) {
281: size_t len;
282: PetscStrlen(dafieldname[s+b],&len);
283: tlen += len + 1; /* field + "-" */
284: }
285: PetscMalloc1(tlen,&fieldname[i]);
286: PetscStrcpy(fieldname[i],"Vector-");
287: for (b=0;b<bss[i]-1;b++) {
288: PetscStrcat(fieldname[i],dafieldname[s+b]);
289: PetscStrcat(fieldname[i],"-");
290: }
291: PetscStrcat(fieldname[i],dafieldname[s+b]);
292: }
293: dims[i] = dim;
294: nlocal[i] = M*N*P*bss[i];
295: s += bss[i];
296: }
298: /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */
299: PetscNew(&ctx);
300: ctx->xlocal = xlocal;
302: /* create work vectors */
303: for (i=0;i<nf;i++) {
304: VecCreateMPI(PetscObjectComm((PetscObject)da),nlocal[i],PETSC_DECIDE,&Ufield[i]);
305: PetscObjectSetName((PetscObject)Ufield[i],fieldname[i]);
306: VecSetBlockSize(Ufield[i],bss[i]);
307: VecSetDM(Ufield[i],da);
308: }
310: PetscViewerGLVisSetFields(viewer,nf,(const char**)fec_type,dims,DMDASampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DMDAFieldDestroyGLVisViewerCtx_Private);
311: for (i=0;i<nf;i++) {
312: PetscFree(fec_type[i]);
313: PetscFree(fieldname[i]);
314: VecDestroy(&Ufield[i]);
315: }
316: PetscFree6(fec_type,nlocal,bss,dims,fieldname,Ufield);
317: }
318: return(0);
319: }
321: static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
322: {
323: DM da,cda;
324: Vec xcoorl;
325: PetscMPIInt size;
326: const PetscScalar *array;
327: PetscContainer glvis_container;
328: PetscInt dim,sdim,i,vid[8],mid,cid,cdof;
329: PetscInt sx,sy,sz,ie,je,ke,ien,jen,ken;
330: PetscInt gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
331: PetscBool enabled = PETSC_TRUE, isascii;
332: PetscErrorCode ierr;
333: const char *fmt;
338: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
339: if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
340: MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);
341: if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
342: DMGetDimension(dm,&dim);
344: /* get container: determines if a process visualizes is portion of the data or not */
345: PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);
346: if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
347: {
348: PetscViewerGLVisInfo glvis_info;
349: PetscContainerGetPointer(glvis_container,(void**)&glvis_info);
350: enabled = glvis_info->enabled;
351: fmt = glvis_info->fmt;
352: }
353: /* this can happen if we are calling DMView outside of VecView_GLVis */
354: PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);
355: if (!da) {DMSetUpGLVisViewer_DMDA((PetscObject)dm,NULL);}
356: PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);
357: if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
358: DMGetCoordinateDim(da,&sdim);
360: PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");
361: PetscViewerASCIIPrintf(viewer,"\ndimension\n");
362: PetscViewerASCIIPrintf(viewer,"%D\n",dim);
364: if (!enabled) {
365: PetscViewerASCIIPrintf(viewer,"\nelements\n");
366: PetscViewerASCIIPrintf(viewer,"%D\n",0);
367: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
368: PetscViewerASCIIPrintf(viewer,"%D\n",0);
369: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
370: PetscViewerASCIIPrintf(viewer,"%D\n",0);
371: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
372: return(0);
373: }
375: DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
376: i = ien;
377: if (dim > 1) i *= jen;
378: if (dim > 2) i *= ken;
379: PetscViewerASCIIPrintf(viewer,"\nelements\n");
380: PetscViewerASCIIPrintf(viewer,"%D\n",i);
381: switch (dim) {
382: case 1:
383: for (ie = 0; ie < ien; ie++) {
384: vid[0] = ie;
385: vid[1] = ie+1;
386: mid = 1; /* material id */
387: cid = 1; /* segment */
388: PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);
389: }
390: break;
391: case 2:
392: for (je = 0; je < jen; je++) {
393: for (ie = 0; ie < ien; ie++) {
394: vid[0] = je*(ien+1) + ie;
395: vid[1] = je*(ien+1) + ie+1;
396: vid[2] = (je+1)*(ien+1) + ie+1;
397: vid[3] = (je+1)*(ien+1) + ie;
398: mid = 1; /* material id */
399: cid = 3; /* quad */
400: PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);
401: }
402: }
403: break;
404: case 3:
405: for (ke = 0; ke < ken; ke++) {
406: for (je = 0; je < jen; je++) {
407: for (ie = 0; ie < ien; ie++) {
408: vid[0] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie;
409: vid[1] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie+1;
410: vid[2] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
411: vid[3] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
412: vid[4] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie;
413: vid[5] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie+1;
414: vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
415: vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
416: mid = 1; /* material id */
417: cid = 5; /* hex */
418: 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]);
419: }
420: }
421: }
422: break;
423: default:
424: SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
425: break;
426: }
427: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
428: PetscViewerASCIIPrintf(viewer,"%D\n",0);
430: /* vertex coordinates */
431: PetscObjectQuery((PetscObject)da,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);
432: if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
433: DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);
434: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
435: PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);
436: VecGetDM(xcoorl,&cda);
437: VecGetArrayRead(xcoorl,&array);
438: if (!cda) { /* HO viz */
439: const char *fecname;
440: PetscInt nc,nl;
442: PetscObjectGetName((PetscObject)xcoorl,&fecname);
443: PetscViewerASCIIPrintf(viewer,"nodes\n");
444: PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");
445: PetscViewerASCIIPrintf(viewer,"FiniteElementCollection: %s\n",fecname);
446: PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);
447: PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n"); /*Ordering::byVDIM*/
448: /* L2 coordinates */
449: DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
450: VecGetLocalSize(xcoorl,&nl);
451: nc = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
452: cdof = nl/nc;
453: if (!ien) ien++;
454: if (!jen) jen++;
455: if (!ken) ken++;
456: ist = jst = kst = 0;
457: gm = ien;
458: gn = jen;
459: gp = ken;
460: } else {
461: DMDAGhostedGLVisViewerCtx *dactx;
463: DMGetApplicationContext(da,(void**)&dactx);
464: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
465: cdof = sdim;
466: DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);
467: DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);
468: if (dactx->ll) {
469: kst = jst = ist = 0;
470: } else {
471: kst = gsz != sz ? 1 : 0;
472: jst = gsy != sy ? 1 : 0;
473: ist = gsx != sx ? 1 : 0;
474: }
475: }
476: for (ke = kst; ke < kst + ken; ke++) {
477: for (je = jst; je < jst + jen; je++) {
478: for (ie = ist; ie < ist + ien; ie++) {
479: PetscInt c;
481: i = ke * gm * gn + je * gm + ie;
482: for (c=0;c<cdof/sdim;c++) {
483: PetscInt d;
484: for (d=0;d<sdim;d++) {
485: PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[cdof*i+c*sdim+d]));
486: }
487: PetscViewerASCIIPrintf(viewer,"\n");
488: }
489: }
490: }
491: }
492: VecRestoreArrayRead(xcoorl,&array);
493: return(0);
494: }
496: /* 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 */
497: PETSC_INTERN PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
498: {
500: PetscBool isglvis,isascii;
505: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);
506: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
507: if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
508: if (isglvis) {
509: PetscViewer view;
510: PetscViewerGLVisType type;
512: PetscViewerGLVisGetType_Private(viewer,&type);
513: PetscViewerGLVisGetDMWindow_Private(viewer,&view);
514: if (view) { /* in the socket case, it may happen that the connection failed */
515: if (type == PETSC_VIEWER_GLVIS_SOCKET) {
516: PetscMPIInt size,rank;
518: MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);
519: MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);
520: PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);
521: }
522: DMDAView_GLVis_ASCII(dm,view);
523: PetscViewerFlush(view);
524: if (type == PETSC_VIEWER_GLVIS_SOCKET) {
525: PetscInt dim;
526: const char* name;
528: PetscObjectGetName((PetscObject)dm,&name);
529: DMGetDimension(dm,&dim);
530: PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);
531: PetscBarrier((PetscObject)dm);
532: }
533: }
534: } else {
535: DMDAView_GLVis_ASCII(dm,viewer);
536: }
537: return(0);
538: }