Actual source code: grglvis.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  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:       break;
198:     }
199:     DMSetApplicationContext(daview,dactx);
200:     DMSetApplicationContextDestroy(daview,DMDAGhostedDestroyGLVisViewerCtx_Private);
201:     DMSetUp(daview);
202:     if (!xcoor) {
203:       DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);
204:       DMGetCoordinates(daview,&xcoor);
205:     }
206:     if (dacoord) {
207:       DMSetUp(dacoord);
208:       DMCreateLocalVector(dacoord,&xcoorl);
209:       DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);
210:       DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);
211:       DMDestroy(&dacoord);
212:     } else {
213:       PetscInt   ien,jen,ken,nc,nl,cdof,deg;
214:       char       fecmesh[64];
215:       const char *name;
216:       PetscBool  flg;

218:       DMDAGetNumElementsGhosted(daview,&ien,&jen,&ken);
219:       nc   = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);

221:       VecGetLocalSize(xcoor,&nl);
222:       if (nc && nl % nc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible local coordinate size %D and number of cells %D",nl,nc);
223:       VecDuplicate(xcoor,&xcoorl);
224:       VecCopy(xcoor,xcoorl);
225:       VecSetDM(xcoorl,NULL);
226:       PetscObjectGetName((PetscObject)xcoor,&name);
227:       PetscStrbeginswith(name,"FiniteElementCollection:",&flg);
228:       if (!flg) {
229:         deg = 0;
230:         if (nc && nl) {
231:           cdof = nl/(nc*dim);
232:           deg  = 1;
233:           while (1) {
234:             PetscInt degd = 1;
235:             for (i=0;i<dim;i++) degd *= (deg+1);
236:             if (degd > cdof) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell dofs %D",cdof);
237:             if (degd == cdof) break;
238:             deg++;
239:           }
240:         }
241:         PetscSNPrintf(fecmesh,sizeof(fecmesh),"FiniteElementCollection: L2_T1_%DD_P%D",dim,deg);
242:         PetscObjectSetName((PetscObject)xcoorl,fecmesh);
243:       } else {
244:         PetscObjectSetName((PetscObject)xcoorl,name);
245:       }
246:     }

248:     /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */
249:     PetscObjectCompose((PetscObject)daview,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);
250:     PetscObjectDereference((PetscObject)xcoorl);

252:     /* daview is composed with the original DMDA */
253:     PetscObjectCompose(oda,"GLVisGraphicsDMDAGhosted",(PetscObject)daview);
254:     PetscObjectDereference((PetscObject)daview);
255:   }

257:   /* customize the viewer if present */
258:   if (viewer) {
259:     DMDAFieldGLVisViewerCtx   *ctx;
260:     DMDAGhostedGLVisViewerCtx *dactx;
261:     char                      fec[64];
262:     Vec                       xlocal,*Ufield;
263:     const char                **dafieldname;
264:     char                      **fec_type,**fieldname;
265:     PetscInt                  *nlocal,*bss,*dims;
266:     PetscInt                  dim,M,N,P,dof,s,i,nf;
267:     PetscBool                 bsset;

269:     DMDAGetInfo(daview,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
270:     DMGetApplicationContext(daview,(void**)&dactx);
271:     DMCreateLocalVector(daview,&xlocal);
272:     DMDAGetFieldNames(da,(const char * const **)&dafieldname);
273:     DMDAGetNumVerticesGhosted(daview,&M,&N,&P);
274:     PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: H1_%DD_P1",dim);
275:     PetscMalloc6(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname,dof,&Ufield);
276:     for (i=0;i<dof;i++) bss[i] = 1;
277:     nf = dof;

279:     PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Field options","PetscViewer");
280:     PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);
281:     PetscOptionsEnd();
282:     if (bsset) {
283:       PetscInt t;
284:       for (i=0,t=0;i<nf;i++) t += bss[i];
285:       if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
286:     } else nf = dof;

288:     for (i=0,s=0;i<nf;i++) {
289:       PetscStrallocpy(fec,&fec_type[i]);
290:       if (bss[i] == 1) {
291:         PetscStrallocpy(dafieldname[s],&fieldname[i]);
292:       } else {
293:         PetscInt b;
294:         size_t tlen = 9; /* "Vector-" + end */
295:         for (b=0;b<bss[i];b++) {
296:           size_t len;
297:           PetscStrlen(dafieldname[s+b],&len);
298:           tlen += len + 1; /* field + "-" */
299:         }
300:         PetscMalloc1(tlen,&fieldname[i]);
301:         PetscStrcpy(fieldname[i],"Vector-");
302:         for (b=0;b<bss[i]-1;b++) {
303:           PetscStrcat(fieldname[i],dafieldname[s+b]);
304:           PetscStrcat(fieldname[i],"-");
305:         }
306:         PetscStrcat(fieldname[i],dafieldname[s+b]);
307:       }
308:       dims[i] = dim;
309:       nlocal[i] = M*N*P*bss[i];
310:       s += bss[i];
311:     }

313:     /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */
314:     PetscNew(&ctx);
315:     ctx->xlocal = xlocal;

317:     /* create work vectors */
318:     for (i=0;i<nf;i++) {
319:       VecCreateMPI(PetscObjectComm((PetscObject)da),nlocal[i],PETSC_DECIDE,&Ufield[i]);
320:       PetscObjectSetName((PetscObject)Ufield[i],fieldname[i]);
321:       VecSetBlockSize(Ufield[i],bss[i]);
322:       VecSetDM(Ufield[i],da);
323:     }

325:     PetscViewerGLVisSetFields(viewer,nf,(const char**)fec_type,dims,DMDASampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DMDAFieldDestroyGLVisViewerCtx_Private);
326:     for (i=0;i<nf;i++) {
327:       PetscFree(fec_type[i]);
328:       PetscFree(fieldname[i]);
329:       VecDestroy(&Ufield[i]);
330:     }
331:     PetscFree6(fec_type,nlocal,bss,dims,fieldname,Ufield);
332:   }
333:   return(0);
334: }

336: static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
337: {
338:   DM                da,cda;
339:   Vec               xcoorl;
340:   PetscMPIInt       size;
341:   const PetscScalar *array;
342:   PetscContainer    glvis_container;
343:   PetscInt          dim,sdim,i,vid[8],mid,cid,cdof;
344:   PetscInt          sx,sy,sz,ie,je,ke,ien,jen,ken,nel;
345:   PetscInt          gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
346:   PetscBool         enabled = PETSC_TRUE, isascii;
347:   PetscErrorCode    ierr;
348:   const char        *fmt;

353:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
354:   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
355:   MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);
356:   if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
357:   DMGetDimension(dm,&dim);

359:   /* get container: determines if a process visualizes is portion of the data or not */
360:   PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);
361:   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
362:   {
363:     PetscViewerGLVisInfo glvis_info;
364:     PetscContainerGetPointer(glvis_container,(void**)&glvis_info);
365:     enabled = glvis_info->enabled;
366:     fmt = glvis_info->fmt;
367:   }
368:   /* this can happen if we are calling DMView outside of VecView_GLVis */
369:   PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);
370:   if (!da) {DMSetUpGLVisViewer_DMDA((PetscObject)dm,NULL);}
371:   PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);
372:   if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
373:   DMGetCoordinateDim(da,&sdim);

375:   PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");
376:   PetscViewerASCIIPrintf(viewer,"\ndimension\n");
377:   PetscViewerASCIIPrintf(viewer,"%D\n",dim);

379:   if (!enabled) {
380:     PetscViewerASCIIPrintf(viewer,"\nelements\n");
381:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
382:     PetscViewerASCIIPrintf(viewer,"\nboundary\n");
383:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
384:     PetscViewerASCIIPrintf(viewer,"\nvertices\n");
385:     PetscViewerASCIIPrintf(viewer,"%D\n",0);
386:     PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
387:     return(0);
388:   }

390:   DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
391:   nel  = ien;
392:   if (dim > 1) nel *= jen;
393:   if (dim > 2) nel *= ken;
394:   PetscViewerASCIIPrintf(viewer,"\nelements\n");
395:   PetscViewerASCIIPrintf(viewer,"%D\n",nel);
396:   switch (dim) {
397:   case 1:
398:     for (ie = 0; ie < ien; ie++) {
399:       vid[0] = ie;
400:       vid[1] = ie+1;
401:       mid    = 1; /* material id */
402:       cid    = 1; /* segment */
403:       PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);
404:     }
405:     break;
406:   case 2:
407:     for (je = 0; je < jen; je++) {
408:       for (ie = 0; ie < ien; ie++) {
409:         vid[0] =     je*(ien+1) + ie;
410:         vid[1] =     je*(ien+1) + ie+1;
411:         vid[2] = (je+1)*(ien+1) + ie+1;
412:         vid[3] = (je+1)*(ien+1) + ie;
413:         mid    = 1; /* material id */
414:         cid    = 3; /* quad */
415:         PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);
416:       }
417:     }
418:     break;
419:   case 3:
420:     for (ke = 0; ke < ken; ke++) {
421:       for (je = 0; je < jen; je++) {
422:         for (ie = 0; ie < ien; ie++) {
423:           vid[0] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie;
424:           vid[1] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
425:           vid[2] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
426:           vid[3] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
427:           vid[4] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie;
428:           vid[5] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
429:           vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
430:           vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
431:           mid    = 1; /* material id */
432:           cid    = 5; /* hex */
433:           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]);
434:         }
435:       }
436:     }
437:     break;
438:   default:
439:     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
440:     break;
441:   }
442:   PetscViewerASCIIPrintf(viewer,"\nboundary\n");
443:   PetscViewerASCIIPrintf(viewer,"%D\n",0);

445:   /* vertex coordinates */
446:   PetscObjectQuery((PetscObject)da,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);
447:   if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
448:   DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);
449:   PetscViewerASCIIPrintf(viewer,"\nvertices\n");
450:   PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);
451:   if (nel) {
452:     VecGetDM(xcoorl,&cda);
453:     VecGetArrayRead(xcoorl,&array);
454:     if (!cda) { /* HO viz */
455:       const char *fecname;
456:       PetscInt   nc,nl;

458:       PetscObjectGetName((PetscObject)xcoorl,&fecname);
459:       PetscViewerASCIIPrintf(viewer,"nodes\n");
460:       PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");
461:       PetscViewerASCIIPrintf(viewer,"%s\n",fecname);
462:       PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);
463:       PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n"); /*Ordering::byVDIM*/
464:       /* L2 coordinates */
465:       DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);
466:       VecGetLocalSize(xcoorl,&nl);
467:       nc   = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
468:       cdof = nc ? nl/nc : 0;
469:       if (!ien) ien++;
470:       if (!jen) jen++;
471:       if (!ken) ken++;
472:       ist = jst = kst = 0;
473:       gm = ien;
474:       gn = jen;
475:       gp = ken;
476:     } else {
477:       DMDAGhostedGLVisViewerCtx *dactx;

479:       DMGetApplicationContext(da,(void**)&dactx);
480:       PetscViewerASCIIPrintf(viewer,"%D\n",sdim);
481:       cdof = sdim;
482:       DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);
483:       DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);
484:       if (dactx->ll) {
485:         kst = jst = ist = 0;
486:       } else {
487:         kst  = gsz != sz ? 1 : 0;
488:         jst  = gsy != sy ? 1 : 0;
489:         ist  = gsx != sx ? 1 : 0;
490:       }
491:     }
492:     for (ke = kst; ke < kst + ken; ke++) {
493:       for (je = jst; je < jst + jen; je++) {
494:         for (ie = ist; ie < ist + ien; ie++) {
495:           PetscInt c;

497:           i = ke * gm * gn + je * gm + ie;
498:           for (c=0;c<cdof/sdim;c++) {
499:             PetscInt d;
500:             for (d=0;d<sdim;d++) {
501:               PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[cdof*i+c*sdim+d]));
502:             }
503:             PetscViewerASCIIPrintf(viewer,"\n");
504:           }
505:         }
506:       }
507:     }
508:     VecRestoreArrayRead(xcoorl,&array);
509:   }
510:   return(0);
511: }

513: PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
514: {
517:   DMView_GLVis(dm,viewer,DMDAView_GLVis_ASCII);
518:   return(0);
519: }