Actual source code: grglvis.c
petsc-3.14.6 2021-03-30
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 = 0,jen = 0,ken = 0,dim;
76: PetscInt tote;
80: DMGetDimension(da,&dim);
81: DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
82: tote = ien * (dim > 1 ? jen : 1) * (dim > 2 ? ken : 1);
83: if (tote) {
84: ien = ien+1;
85: jen = dim > 1 ? jen+1 : 1;
86: ken = dim > 2 ? ken+1 : 1;
87: }
88: if (nvx) *nvx = ien;
89: if (nvy) *nvy = jen;
90: if (nvz) *nvz = ken;
91: return(0);
92: }
94: static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx)
95: {
96: DM da;
97: DMDAFieldGLVisViewerCtx *ctx = (DMDAFieldGLVisViewerCtx*)vctx;
98: DMDAGhostedGLVisViewerCtx *dactx;
99: const PetscScalar *array;
100: PetscScalar **arrayf;
101: PetscInt i,f,ii,ien,jen,ken,ie,je,ke,bs,*bss;
102: PetscInt sx,sy,sz,gsx,gsy,gsz,ist,jst,kst,gm,gn,gp;
103: PetscErrorCode ierr;
106: VecGetDM(ctx->xlocal,&da);
107: if (!da) SETERRQ(PetscObjectComm(oX),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
108: DMGetApplicationContext(da,(void**)&dactx);
109: VecGetBlockSize(ctx->xlocal,&bs);
110: DMGlobalToLocalBegin(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);
111: DMGlobalToLocalEnd(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);
112: DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);
113: DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);
114: DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);
115: if (dactx->ll) {
116: kst = jst = ist = 0;
117: } else {
118: kst = gsz != sz ? 1 : 0;
119: jst = gsy != sy ? 1 : 0;
120: ist = gsx != sx ? 1 : 0;
121: }
122: PetscMalloc2(nf,&arrayf,nf,&bss);
123: VecGetArrayRead(ctx->xlocal,&array);
124: for (f=0;f<nf;f++) {
125: VecGetBlockSize((Vec)oXf[f],&bss[f]);
126: VecGetArray((Vec)oXf[f],&arrayf[f]);
127: }
128: for (ke = kst, ii = 0; ke < kst + ken; ke++) {
129: for (je = jst; je < jst + jen; je++) {
130: for (ie = ist; ie < ist + ien; ie++) {
131: PetscInt cf,b;
132: i = ke * gm * gn + je * gm + ie;
133: for (f=0,cf=0;f<nf;f++)
134: for (b=0;b<bss[f];b++)
135: arrayf[f][bss[f]*ii+b] = array[i*bs+cf++];
136: ii++;
137: }
138: }
139: }
140: for (f=0;f<nf;f++) { VecRestoreArray((Vec)oXf[f],&arrayf[f]); }
141: VecRestoreArrayRead(ctx->xlocal,&array);
142: PetscFree2(arrayf,bss);
143: return(0);
144: }
146: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
147: {
148: DM da = (DM)oda,daview;
152: PetscObjectQuery(oda,"GLVisGraphicsDMDAGhosted",(PetscObject*)&daview);
153: if (!daview) {
154: DMDAGhostedGLVisViewerCtx *dactx;
155: DM dacoord = NULL;
156: Vec xcoor,xcoorl;
157: PetscBool hashocoord = PETSC_FALSE;
158: const PetscInt *lx,*ly,*lz;
159: PetscInt dim,M,N,P,m,n,p,dof,s,i;
161: PetscNew(&dactx);
162: dactx->ll = PETSC_TRUE; /* default to match elements layout obtained by DMDAGetElements */
163: PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");
164: PetscOptionsBool("-viewer_glvis_dm_da_ll","Left-looking subdomain view",NULL,dactx->ll,&dactx->ll,NULL);
165: PetscOptionsEnd();
166: /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
167: DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL);
168: DMDAGetOwnershipRanges(da,&lx,&ly,&lz);
169: PetscObjectQuery((PetscObject)da,"_glvis_mesh_coords",(PetscObject*)&xcoor);
170: if (!xcoor) {
171: DMGetCoordinates(da,&xcoor);
172: } else {
173: hashocoord = PETSC_TRUE;
174: }
175: PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n");
176: switch (dim) {
177: case 1:
178: DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview);
179: if (!hashocoord) {
180: DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord);
181: }
182: break;
183: case 2:
184: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,dof,1,lx,ly,&daview);
185: if (!hashocoord) {
186: DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,2,1,lx,ly,&dacoord);
187: }
188: break;
189: case 3:
190: 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);
191: if (!hashocoord) {
192: 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);
193: }
194: break;
195: default:
196: SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
197: }
198: DMSetApplicationContext(daview,dactx);
199: DMSetApplicationContextDestroy(daview,DMDAGhostedDestroyGLVisViewerCtx_Private);
200: DMSetUp(daview);
201: if (!xcoor) {
202: DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);
203: DMGetCoordinates(daview,&xcoor);
204: }
205: if (dacoord) {
206: DMSetUp(dacoord);
207: DMCreateLocalVector(dacoord,&xcoorl);
208: DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);
209: DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);
210: DMDestroy(&dacoord);
211: } else {
212: PetscInt ien,jen,ken,nc,nl,cdof,deg;
213: char fecmesh[64];
214: const char *name;
215: PetscBool flg;
217: DMDAGetNumElementsGhosted(daview,&ien,&jen,&ken);
218: nc = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
220: VecGetLocalSize(xcoor,&nl);
221: if (nc && nl % nc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible local coordinate size %D and number of cells %D",nl,nc);
222: VecDuplicate(xcoor,&xcoorl);
223: VecCopy(xcoor,xcoorl);
224: VecSetDM(xcoorl,NULL);
225: PetscObjectGetName((PetscObject)xcoor,&name);
226: PetscStrbeginswith(name,"FiniteElementCollection:",&flg);
227: if (!flg) {
228: deg = 0;
229: if (nc && nl) {
230: cdof = nl/(nc*dim);
231: deg = 1;
232: while (1) {
233: PetscInt degd = 1;
234: for (i=0;i<dim;i++) degd *= (deg+1);
235: if (degd > cdof) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell dofs %D",cdof);
236: if (degd == cdof) break;
237: deg++;
238: }
239: }
240: PetscSNPrintf(fecmesh,sizeof(fecmesh),"FiniteElementCollection: L2_T1_%DD_P%D",dim,deg);
241: PetscObjectSetName((PetscObject)xcoorl,fecmesh);
242: } else {
243: PetscObjectSetName((PetscObject)xcoorl,name);
244: }
245: }
247: /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */
248: PetscObjectCompose((PetscObject)daview,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);
249: PetscObjectDereference((PetscObject)xcoorl);
251: /* daview is composed with the original DMDA */
252: PetscObjectCompose(oda,"GLVisGraphicsDMDAGhosted",(PetscObject)daview);
253: PetscObjectDereference((PetscObject)daview);
254: }
256: /* customize the viewer if present */
257: if (viewer) {
258: DMDAFieldGLVisViewerCtx *ctx;
259: DMDAGhostedGLVisViewerCtx *dactx;
260: char fec[64];
261: Vec xlocal,*Ufield;
262: const char **dafieldname;
263: char **fec_type,**fieldname;
264: PetscInt *nlocal,*bss,*dims;
265: PetscInt dim,M,N,P,dof,s,i,nf;
266: PetscBool bsset;
268: DMDAGetInfo(daview,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
269: DMGetApplicationContext(daview,(void**)&dactx);
270: DMCreateLocalVector(daview,&xlocal);
271: DMDAGetFieldNames(da,(const char * const **)&dafieldname);
272: DMDAGetNumVerticesGhosted(daview,&M,&N,&P);
273: PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: H1_%DD_P1",dim);
274: PetscMalloc6(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname,dof,&Ufield);
275: for (i=0;i<dof;i++) bss[i] = 1;
276: nf = dof;
278: PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Field options","PetscViewer");
279: PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);
280: PetscOptionsEnd();
281: if (bsset) {
282: PetscInt t;
283: for (i=0,t=0;i<nf;i++) t += bss[i];
284: if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
285: } else nf = dof;
287: for (i=0,s=0;i<nf;i++) {
288: PetscStrallocpy(fec,&fec_type[i]);
289: if (bss[i] == 1) {
290: PetscStrallocpy(dafieldname[s],&fieldname[i]);
291: } else {
292: PetscInt b;
293: size_t tlen = 9; /* "Vector-" + end */
294: for (b=0;b<bss[i];b++) {
295: size_t len;
296: PetscStrlen(dafieldname[s+b],&len);
297: tlen += len + 1; /* field + "-" */
298: }
299: PetscMalloc1(tlen,&fieldname[i]);
300: PetscStrcpy(fieldname[i],"Vector-");
301: for (b=0;b<bss[i]-1;b++) {
302: PetscStrcat(fieldname[i],dafieldname[s+b]);
303: PetscStrcat(fieldname[i],"-");
304: }
305: PetscStrcat(fieldname[i],dafieldname[s+b]);
306: }
307: dims[i] = dim;
308: nlocal[i] = M*N*P*bss[i];
309: s += bss[i];
310: }
312: /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */
313: PetscNew(&ctx);
314: ctx->xlocal = xlocal;
316: /* create work vectors */
317: for (i=0;i<nf;i++) {
318: VecCreateMPI(PetscObjectComm((PetscObject)da),nlocal[i],PETSC_DECIDE,&Ufield[i]);
319: PetscObjectSetName((PetscObject)Ufield[i],fieldname[i]);
320: VecSetBlockSize(Ufield[i],bss[i]);
321: VecSetDM(Ufield[i],da);
322: }
324: PetscViewerGLVisSetFields(viewer,nf,(const char**)fec_type,dims,DMDASampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DMDAFieldDestroyGLVisViewerCtx_Private);
325: for (i=0;i<nf;i++) {
326: PetscFree(fec_type[i]);
327: PetscFree(fieldname[i]);
328: VecDestroy(&Ufield[i]);
329: }
330: PetscFree6(fec_type,nlocal,bss,dims,fieldname,Ufield);
331: }
332: return(0);
333: }
335: static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
336: {
337: DM da,cda;
338: Vec xcoorl;
339: PetscMPIInt size;
340: const PetscScalar *array;
341: PetscContainer glvis_container;
342: PetscInt dim,sdim,i,vid[8],mid,cid,cdof;
343: PetscInt sx,sy,sz,ie,je,ke,ien,jen,ken,nel;
344: PetscInt gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
345: PetscBool enabled = PETSC_TRUE, isascii;
346: PetscErrorCode ierr;
347: const char *fmt;
352: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
353: if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
354: MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);
355: if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
356: DMGetDimension(dm,&dim);
358: /* get container: determines if a process visualizes is portion of the data or not */
359: PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);
360: if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
361: {
362: PetscViewerGLVisInfo glvis_info;
363: PetscContainerGetPointer(glvis_container,(void**)&glvis_info);
364: enabled = glvis_info->enabled;
365: fmt = glvis_info->fmt;
366: }
367: /* this can happen if we are calling DMView outside of VecView_GLVis */
368: PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);
369: if (!da) {DMSetUpGLVisViewer_DMDA((PetscObject)dm,NULL);}
370: PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);
371: if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
372: DMGetCoordinateDim(da,&sdim);
374: PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");
375: PetscViewerASCIIPrintf(viewer,"\ndimension\n");
376: PetscViewerASCIIPrintf(viewer,"%D\n",dim);
378: if (!enabled) {
379: PetscViewerASCIIPrintf(viewer,"\nelements\n");
380: PetscViewerASCIIPrintf(viewer,"%D\n",0);
381: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
382: PetscViewerASCIIPrintf(viewer,"%D\n",0);
383: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
384: PetscViewerASCIIPrintf(viewer,"%D\n",0);
385: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
386: return(0);
387: }
389: DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
390: nel = ien;
391: if (dim > 1) nel *= jen;
392: if (dim > 2) nel *= ken;
393: PetscViewerASCIIPrintf(viewer,"\nelements\n");
394: PetscViewerASCIIPrintf(viewer,"%D\n",nel);
395: switch (dim) {
396: case 1:
397: for (ie = 0; ie < ien; ie++) {
398: vid[0] = ie;
399: vid[1] = ie+1;
400: mid = 1; /* material id */
401: cid = 1; /* segment */
402: PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);
403: }
404: break;
405: case 2:
406: for (je = 0; je < jen; je++) {
407: for (ie = 0; ie < ien; ie++) {
408: vid[0] = je*(ien+1) + ie;
409: vid[1] = je*(ien+1) + ie+1;
410: vid[2] = (je+1)*(ien+1) + ie+1;
411: vid[3] = (je+1)*(ien+1) + ie;
412: mid = 1; /* material id */
413: cid = 3; /* quad */
414: PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);
415: }
416: }
417: break;
418: case 3:
419: for (ke = 0; ke < ken; ke++) {
420: for (je = 0; je < jen; je++) {
421: for (ie = 0; ie < ien; ie++) {
422: vid[0] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie;
423: vid[1] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie+1;
424: vid[2] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
425: vid[3] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
426: vid[4] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie;
427: vid[5] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie+1;
428: vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
429: vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
430: mid = 1; /* material id */
431: cid = 5; /* hex */
432: 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]);
433: }
434: }
435: }
436: break;
437: default:
438: SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
439: }
440: PetscViewerASCIIPrintf(viewer,"\nboundary\n");
441: PetscViewerASCIIPrintf(viewer,"%D\n",0);
443: /* vertex coordinates */
444: PetscObjectQuery((PetscObject)da,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);
445: if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
446: DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);
447: PetscViewerASCIIPrintf(viewer,"\nvertices\n");
448: PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);
449: if (nel) {
450: VecGetDM(xcoorl,&cda);
451: VecGetArrayRead(xcoorl,&array);
452: if (!cda) { /* HO viz */
453: const char *fecname;
454: PetscInt nc,nl;
456: PetscObjectGetName((PetscObject)xcoorl,&fecname);
457: PetscViewerASCIIPrintf(viewer,"nodes\n");
458: PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");
459: PetscViewerASCIIPrintf(viewer,"%s\n",fecname);
460: PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);
461: PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n"); /*Ordering::byVDIM*/
462: /* L2 coordinates */
463: DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
464: VecGetLocalSize(xcoorl,&nl);
465: nc = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
466: cdof = nc ? nl/nc : 0;
467: if (!ien) ien++;
468: if (!jen) jen++;
469: if (!ken) ken++;
470: ist = jst = kst = 0;
471: gm = ien;
472: gn = jen;
473: gp = ken;
474: } else {
475: DMDAGhostedGLVisViewerCtx *dactx;
477: DMGetApplicationContext(da,(void**)&dactx);
478: PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
479: cdof = sdim;
480: DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);
481: DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);
482: if (dactx->ll) {
483: kst = jst = ist = 0;
484: } else {
485: kst = gsz != sz ? 1 : 0;
486: jst = gsy != sy ? 1 : 0;
487: ist = gsx != sx ? 1 : 0;
488: }
489: }
490: for (ke = kst; ke < kst + ken; ke++) {
491: for (je = jst; je < jst + jen; je++) {
492: for (ie = ist; ie < ist + ien; ie++) {
493: PetscInt c;
495: i = ke * gm * gn + je * gm + ie;
496: for (c=0;c<cdof/sdim;c++) {
497: PetscInt d;
498: for (d=0;d<sdim;d++) {
499: PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[cdof*i+c*sdim+d]));
500: }
501: PetscViewerASCIIPrintf(viewer,"\n");
502: }
503: }
504: }
505: }
506: VecRestoreArrayRead(xcoorl,&array);
507: }
508: return(0);
509: }
511: PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
512: {
515: DMView_GLVis(dm,viewer,DMDAView_GLVis_ASCII);
516: return(0);
517: }