xref: /petsc/src/dm/impls/plex/plexglvis.c (revision 8135c375e96a7bf70d67c02f02cd826154817056)
1*8135c375SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
2*8135c375SStefano Zampini #include <petsc/private/petscimpl.h>
3*8135c375SStefano Zampini #include <petsc/private/dmpleximpl.h>
4*8135c375SStefano Zampini #include <petscbt.h>
5*8135c375SStefano Zampini #include <petscdmplex.h>
6*8135c375SStefano Zampini #include <petscsf.h>
7*8135c375SStefano Zampini #include <petscds.h>
8*8135c375SStefano Zampini 
9*8135c375SStefano Zampini typedef struct {
10*8135c375SStefano Zampini   PetscInt   nf;
11*8135c375SStefano Zampini   VecScatter *scctx;
12*8135c375SStefano Zampini } GLVisViewerCtx;
13*8135c375SStefano Zampini 
14*8135c375SStefano Zampini static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx)
15*8135c375SStefano Zampini {
16*8135c375SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
17*8135c375SStefano Zampini   PetscInt       i;
18*8135c375SStefano Zampini   PetscErrorCode ierr;
19*8135c375SStefano Zampini 
20*8135c375SStefano Zampini   PetscFunctionBegin;
21*8135c375SStefano Zampini   for (i=0;i<ctx->nf;i++) {
22*8135c375SStefano Zampini     ierr = VecScatterDestroy(&ctx->scctx[i]);CHKERRQ(ierr);
23*8135c375SStefano Zampini   }
24*8135c375SStefano Zampini   ierr = PetscFree(ctx->scctx);CHKERRQ(ierr);
25*8135c375SStefano Zampini   ierr = PetscFree(vctx);
26*8135c375SStefano Zampini   PetscFunctionReturn(0);
27*8135c375SStefano Zampini }
28*8135c375SStefano Zampini 
29*8135c375SStefano Zampini static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx)
30*8135c375SStefano Zampini {
31*8135c375SStefano Zampini   GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx;
32*8135c375SStefano Zampini   PetscInt       f;
33*8135c375SStefano Zampini   PetscErrorCode ierr;
34*8135c375SStefano Zampini 
35*8135c375SStefano Zampini   PetscFunctionBegin;
36*8135c375SStefano Zampini   for (f=0;f<nf;f++) {
37*8135c375SStefano Zampini     ierr = VecScatterBegin(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
38*8135c375SStefano Zampini     ierr = VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
39*8135c375SStefano Zampini   }
40*8135c375SStefano Zampini   PetscFunctionReturn(0);
41*8135c375SStefano Zampini }
42*8135c375SStefano Zampini 
43*8135c375SStefano Zampini /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */
44*8135c375SStefano Zampini PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer)
45*8135c375SStefano Zampini {
46*8135c375SStefano Zampini   DM             dm = (DM)odm;
47*8135c375SStefano Zampini   Vec            xlocal,xfield;
48*8135c375SStefano Zampini   PetscDS        ds;
49*8135c375SStefano Zampini   IS             globalNum,isfield;
50*8135c375SStefano Zampini   PetscBT        vown;
51*8135c375SStefano Zampini   char           **fieldname = NULL,**fec_type = NULL;
52*8135c375SStefano Zampini   const PetscInt *gNum;
53*8135c375SStefano Zampini   PetscInt       *nlocal,*bs,*idxs;
54*8135c375SStefano Zampini   PetscInt       f,maxfields,nfields,c,totc,totdofs,Nv,cum,i;
55*8135c375SStefano Zampini   PetscInt       dim,depth,cStart,cEnd,cEndInterior,vStart,vEnd,coneSize;
56*8135c375SStefano Zampini   GLVisViewerCtx *ctx;
57*8135c375SStefano Zampini   PetscSection   s;
58*8135c375SStefano Zampini   PetscBool      simplex;
59*8135c375SStefano Zampini   PetscErrorCode ierr;
60*8135c375SStefano Zampini 
61*8135c375SStefano Zampini   PetscFunctionBegin;
62*8135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
63*8135c375SStefano Zampini   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
64*8135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
65*8135c375SStefano Zampini   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
66*8135c375SStefano Zampini   ierr = DMPlexGetConeSize(dm,cStart,&coneSize);CHKERRQ(ierr);
67*8135c375SStefano Zampini   if (depth == 1) {
68*8135c375SStefano Zampini     if      (coneSize == dim+1)    simplex = PETSC_TRUE;
69*8135c375SStefano Zampini     else if (coneSize == 1 << dim) simplex = PETSC_FALSE;
70*8135c375SStefano Zampini     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells");
71*8135c375SStefano Zampini   }
72*8135c375SStefano Zampini   else if (depth == dim) {
73*8135c375SStefano Zampini     if      (coneSize == dim+1)   simplex = PETSC_TRUE;
74*8135c375SStefano Zampini     else if (coneSize == 2 * dim) simplex = PETSC_FALSE;
75*8135c375SStefano Zampini     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells");
76*8135c375SStefano Zampini   }
77*8135c375SStefano Zampini   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support cell-vertex meshes or interpolated meshes");
78*8135c375SStefano Zampini 
79*8135c375SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
80*8135c375SStefano Zampini   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
81*8135c375SStefano Zampini   ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
82*8135c375SStefano Zampini   ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
83*8135c375SStefano Zampini   ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr);
84*8135c375SStefano Zampini   for (c = cStart, totc = 0; c < cEnd; c++) {
85*8135c375SStefano Zampini     if (gNum[c-cStart] >= 0) {
86*8135c375SStefano Zampini       PetscInt i,numPoints,*points = NULL;
87*8135c375SStefano Zampini 
88*8135c375SStefano Zampini       totc++;
89*8135c375SStefano Zampini       ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
90*8135c375SStefano Zampini       for (i=0;i<numPoints*2;i+= 2) {
91*8135c375SStefano Zampini         if ((points[i] >= vStart) && (points[i] < vEnd)) {
92*8135c375SStefano Zampini           ierr = PetscBTSet(vown,points[i]-vStart);CHKERRQ(ierr);
93*8135c375SStefano Zampini         }
94*8135c375SStefano Zampini       }
95*8135c375SStefano Zampini       ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
96*8135c375SStefano Zampini     }
97*8135c375SStefano Zampini   }
98*8135c375SStefano Zampini   for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscBTLookup(vown,f)) Nv++;
99*8135c375SStefano Zampini 
100*8135c375SStefano Zampini   ierr = DMCreateLocalVector(dm,&xlocal);CHKERRQ(ierr);
101*8135c375SStefano Zampini   ierr = VecGetLocalSize(xlocal,&totdofs);CHKERRQ(ierr);
102*8135c375SStefano Zampini   ierr = DMGetDefaultSection(dm,&s);CHKERRQ(ierr);
103*8135c375SStefano Zampini   ierr = PetscSectionGetNumFields(s,&nfields);CHKERRQ(ierr);
104*8135c375SStefano Zampini   for (f=0,maxfields=0;f<nfields;f++) {
105*8135c375SStefano Zampini     PetscInt bs;
106*8135c375SStefano Zampini 
107*8135c375SStefano Zampini     ierr = PetscSectionGetFieldComponents(s,f,&bs);CHKERRQ(ierr);
108*8135c375SStefano Zampini     maxfields += bs;
109*8135c375SStefano Zampini   }
110*8135c375SStefano Zampini   ierr = PetscCalloc5(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&fec_type,totdofs,&idxs);CHKERRQ(ierr);
111*8135c375SStefano Zampini   ierr = PetscNew(&ctx);CHKERRQ(ierr);
112*8135c375SStefano Zampini   ierr = PetscCalloc1(maxfields,&ctx->scctx);CHKERRQ(ierr);
113*8135c375SStefano Zampini   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
114*8135c375SStefano Zampini   if (ds) {
115*8135c375SStefano Zampini     for (f=0;f<nfields;f++) {
116*8135c375SStefano Zampini       const char* fname;
117*8135c375SStefano Zampini       char        name[256];
118*8135c375SStefano Zampini       PetscObject disc;
119*8135c375SStefano Zampini       size_t      len;
120*8135c375SStefano Zampini 
121*8135c375SStefano Zampini       ierr = PetscSectionGetFieldName(s,f,&fname);CHKERRQ(ierr);
122*8135c375SStefano Zampini       ierr = PetscStrlen(fname,&len);CHKERRQ(ierr);
123*8135c375SStefano Zampini       if (len) {
124*8135c375SStefano Zampini         ierr = PetscStrcpy(name,fname);CHKERRQ(ierr);
125*8135c375SStefano Zampini       } else {
126*8135c375SStefano Zampini         ierr = PetscSNPrintf(name,256,"Field%d",f);CHKERRQ(ierr);
127*8135c375SStefano Zampini       }
128*8135c375SStefano Zampini       ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr);
129*8135c375SStefano Zampini       if (disc) {
130*8135c375SStefano Zampini         PetscClassId id;
131*8135c375SStefano Zampini         PetscInt     Nc;
132*8135c375SStefano Zampini         char         fec[64];
133*8135c375SStefano Zampini 
134*8135c375SStefano Zampini         ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
135*8135c375SStefano Zampini         if (id == PETSCFE_CLASSID) {
136*8135c375SStefano Zampini           PetscFE            fem = (PetscFE)disc;
137*8135c375SStefano Zampini           PetscDualSpace     sp;
138*8135c375SStefano Zampini           PetscDualSpaceType spname;
139*8135c375SStefano Zampini           PetscInt           order;
140*8135c375SStefano Zampini           PetscBool          islag,continuous,H1 = PETSC_TRUE;
141*8135c375SStefano Zampini 
142*8135c375SStefano Zampini           ierr = PetscFEGetNumComponents(fem,&Nc);CHKERRQ(ierr);
143*8135c375SStefano Zampini           ierr = PetscFEGetDualSpace(fem,&sp);CHKERRQ(ierr);
144*8135c375SStefano Zampini           ierr = PetscDualSpaceGetType(sp,&spname);CHKERRQ(ierr);
145*8135c375SStefano Zampini           ierr = PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);CHKERRQ(ierr);
146*8135c375SStefano Zampini           if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space");
147*8135c375SStefano Zampini           ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr);
148*8135c375SStefano Zampini           if (!continuous) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported");
149*8135c375SStefano Zampini           ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr);
150*8135c375SStefano Zampini           if (continuous && order > 0) {
151*8135c375SStefano Zampini             ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%dD_P1",dim);CHKERRQ(ierr);
152*8135c375SStefano Zampini           } else {
153*8135c375SStefano Zampini             H1   = PETSC_FALSE;
154*8135c375SStefano Zampini             ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P%d",dim,order);CHKERRQ(ierr);
155*8135c375SStefano Zampini           }
156*8135c375SStefano Zampini           ierr = PetscStrallocpy(name,&fieldname[ctx->nf]);CHKERRQ(ierr);
157*8135c375SStefano Zampini           bs[ctx->nf] = Nc;
158*8135c375SStefano Zampini           if (H1) {
159*8135c375SStefano Zampini             nlocal[ctx->nf] = Nc * Nv;
160*8135c375SStefano Zampini             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
161*8135c375SStefano Zampini             ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr);
162*8135c375SStefano Zampini             for (i=0,cum=0;i<vEnd-vStart;i++) {
163*8135c375SStefano Zampini               PetscInt j,off;
164*8135c375SStefano Zampini 
165*8135c375SStefano Zampini               if (PetscUnlikely(!PetscBTLookup(vown,i))) continue;
166*8135c375SStefano Zampini               ierr = PetscSectionGetFieldOffset(s,i+vStart,f,&off);CHKERRQ(ierr);
167*8135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
168*8135c375SStefano Zampini             }
169*8135c375SStefano Zampini             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
170*8135c375SStefano Zampini           } else {
171*8135c375SStefano Zampini             nlocal[ctx->nf] = Nc * totc;
172*8135c375SStefano Zampini             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
173*8135c375SStefano Zampini             ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr);
174*8135c375SStefano Zampini             for (i=0,cum=0;i<cEnd-cStart;i++) {
175*8135c375SStefano Zampini               PetscInt j,off;
176*8135c375SStefano Zampini 
177*8135c375SStefano Zampini               if (PetscUnlikely(gNum[i] < 0)) continue;
178*8135c375SStefano Zampini               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
179*8135c375SStefano Zampini               for (j=0;j<Nc;j++) idxs[cum++] = off + j;
180*8135c375SStefano Zampini             }
181*8135c375SStefano Zampini             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
182*8135c375SStefano Zampini           }
183*8135c375SStefano Zampini           ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
184*8135c375SStefano Zampini           ierr = VecDestroy(&xfield);CHKERRQ(ierr);
185*8135c375SStefano Zampini           ierr = ISDestroy(&isfield);CHKERRQ(ierr);
186*8135c375SStefano Zampini           ctx->nf++;
187*8135c375SStefano Zampini         } else if (id == PETSCFV_CLASSID) {
188*8135c375SStefano Zampini           PetscInt c;
189*8135c375SStefano Zampini 
190*8135c375SStefano Zampini           ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr);
191*8135c375SStefano Zampini           ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P0",dim);CHKERRQ(ierr);
192*8135c375SStefano Zampini           for (c = 0; c < Nc; c++) {
193*8135c375SStefano Zampini             char comp[256];
194*8135c375SStefano Zampini             ierr = PetscSNPrintf(comp,256,"%s-Comp%d",name,c);CHKERRQ(ierr);
195*8135c375SStefano Zampini             ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr);
196*8135c375SStefano Zampini             bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */
197*8135c375SStefano Zampini             nlocal[ctx->nf] = totc;
198*8135c375SStefano Zampini             ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr);
199*8135c375SStefano Zampini             ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr);
200*8135c375SStefano Zampini             for (i=0,cum=0;i<cEnd-cStart;i++) {
201*8135c375SStefano Zampini               PetscInt off;
202*8135c375SStefano Zampini 
203*8135c375SStefano Zampini               if (PetscUnlikely(gNum[i])<0) continue;
204*8135c375SStefano Zampini               ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr);
205*8135c375SStefano Zampini               idxs[cum++] = off + c;
206*8135c375SStefano Zampini             }
207*8135c375SStefano Zampini             ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr);
208*8135c375SStefano Zampini             ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr);
209*8135c375SStefano Zampini             ierr = VecDestroy(&xfield);CHKERRQ(ierr);
210*8135c375SStefano Zampini             ierr = ISDestroy(&isfield);CHKERRQ(ierr);
211*8135c375SStefano Zampini             ctx->nf++;
212*8135c375SStefano Zampini           }
213*8135c375SStefano Zampini         } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f);
214*8135c375SStefano Zampini       } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f);
215*8135c375SStefano Zampini     }
216*8135c375SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM");
217*8135c375SStefano Zampini   ierr = PetscBTDestroy(&vown);CHKERRQ(ierr);
218*8135c375SStefano Zampini   ierr = VecDestroy(&xlocal);CHKERRQ(ierr);
219*8135c375SStefano Zampini   ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
220*8135c375SStefano Zampini 
221*8135c375SStefano Zampini   /* customize the viewer */
222*8135c375SStefano Zampini   ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fieldname,(const char**)fec_type,nlocal,bs,DMPlexSampleGLVisFields_Private,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
223*8135c375SStefano Zampini   for (f=0;f<ctx->nf;f++) {
224*8135c375SStefano Zampini     ierr = PetscFree(fieldname[f]);CHKERRQ(ierr);
225*8135c375SStefano Zampini     ierr = PetscFree(fec_type[f]);CHKERRQ(ierr);
226*8135c375SStefano Zampini   }
227*8135c375SStefano Zampini   ierr = PetscFree5(fieldname,nlocal,bs,fec_type,idxs);CHKERRQ(ierr);
228*8135c375SStefano Zampini   PetscFunctionReturn(0);
229*8135c375SStefano Zampini }
230*8135c375SStefano Zampini 
231*8135c375SStefano Zampini typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_UNDEF} MFEM_cid;
232*8135c375SStefano Zampini 
233*8135c375SStefano Zampini MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF  ,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF, MFEM_UNDEF},
234*8135c375SStefano Zampini                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_UNDEF      ,MFEM_UNDEF, MFEM_UNDEF},
235*8135c375SStefano Zampini                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE     ,MFEM_UNDEF, MFEM_UNDEF},
236*8135c375SStefano Zampini                                   {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF   ,MFEM_TETRAHEDRON,MFEM_UNDEF, MFEM_CUBE } };
237*8135c375SStefano Zampini 
238*8135c375SStefano Zampini #undef __FUNCT__
239*8135c375SStefano Zampini #define __FUNCT__ "DMPlexGetPointMFEMCellID_Internal"
240*8135c375SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt p, PetscInt *mid, PetscInt *cid)
241*8135c375SStefano Zampini {
242*8135c375SStefano Zampini   DMLabel        dlabel;
243*8135c375SStefano Zampini   PetscInt       depth,csize;
244*8135c375SStefano Zampini   PetscErrorCode ierr;
245*8135c375SStefano Zampini 
246*8135c375SStefano Zampini   PetscFunctionBegin;
247*8135c375SStefano Zampini   ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
248*8135c375SStefano Zampini   ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
249*8135c375SStefano Zampini   ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr);
250*8135c375SStefano Zampini   if (label) {
251*8135c375SStefano Zampini     ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr);
252*8135c375SStefano Zampini   } else {
253*8135c375SStefano Zampini     *mid = 1;
254*8135c375SStefano Zampini   }
255*8135c375SStefano Zampini   *cid = mfem_table_cid[depth][csize];
256*8135c375SStefano Zampini   PetscFunctionReturn(0);
257*8135c375SStefano Zampini }
258*8135c375SStefano Zampini 
259*8135c375SStefano Zampini #undef __FUNCT__
260*8135c375SStefano Zampini #define __FUNCT__ "DMPlexGetPointMFEMVertexIDs_Internal"
261*8135c375SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscInt *nv, PetscInt vids[])
262*8135c375SStefano Zampini {
263*8135c375SStefano Zampini   PetscInt       dim,i,q,vStart,vEnd,numPoints,*points = NULL;
264*8135c375SStefano Zampini   PetscErrorCode ierr;
265*8135c375SStefano Zampini 
266*8135c375SStefano Zampini   PetscFunctionBegin;
267*8135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
268*8135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
269*8135c375SStefano Zampini   ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
270*8135c375SStefano Zampini   for (i=0,q=0;i<numPoints*2;i+= 2)
271*8135c375SStefano Zampini     if ((points[i] >= vStart) && (points[i] < vEnd))
272*8135c375SStefano Zampini       vids[q++] = points[i]-vStart;
273*8135c375SStefano Zampini   ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr);
274*8135c375SStefano Zampini   ierr = DMPlexInvertCell(dim,q,vids);CHKERRQ(ierr);
275*8135c375SStefano Zampini   *nv  = q;
276*8135c375SStefano Zampini   PetscFunctionReturn(0);
277*8135c375SStefano Zampini }
278*8135c375SStefano Zampini 
279*8135c375SStefano Zampini /* ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR */
280*8135c375SStefano Zampini static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer)
281*8135c375SStefano Zampini {
282*8135c375SStefano Zampini   DM                   cdm;
283*8135c375SStefano Zampini   DMLabel              label;
284*8135c375SStefano Zampini   PetscSection         coordSection,parentSection;
285*8135c375SStefano Zampini   Vec                  coordinates;
286*8135c375SStefano Zampini   IS                   globalNum = NULL;
287*8135c375SStefano Zampini   const PetscScalar    *array;
288*8135c375SStefano Zampini   const PetscInt       *gNum;
289*8135c375SStefano Zampini   PetscInt             bf,p,sdim,dim,depth,novl;
290*8135c375SStefano Zampini   PetscInt             cStart,cEnd,cEndInterior,fStart,fEnd,fEndInterior,vStart,vEnd;
291*8135c375SStefano Zampini   PetscMPIInt          commsize;
292*8135c375SStefano Zampini   PetscBool            ovl,isascii,enable_boundary,enable_ncmesh;
293*8135c375SStefano Zampini   PetscBT              pown;
294*8135c375SStefano Zampini   PetscErrorCode       ierr;
295*8135c375SStefano Zampini   PetscInt             *faces,fpc,vpf;
296*8135c375SStefano Zampini   PetscInt             fv1[]     = {0,1},
297*8135c375SStefano Zampini                        fv2tri[]  = {0,1,
298*8135c375SStefano Zampini                                     1,2,
299*8135c375SStefano Zampini                                     2,0},
300*8135c375SStefano Zampini                        fv2quad[] = {0,1,
301*8135c375SStefano Zampini                                     1,2,
302*8135c375SStefano Zampini                                     2,3,
303*8135c375SStefano Zampini                                     3,0},
304*8135c375SStefano Zampini                        fv3tet[]  = {0,1,2,
305*8135c375SStefano Zampini                                     0,3,1,
306*8135c375SStefano Zampini                                     0,2,3,
307*8135c375SStefano Zampini                                     2,1,3},
308*8135c375SStefano Zampini                        fv3hex[]  = {0,1,2,3,
309*8135c375SStefano Zampini                                  4,5,6,7,
310*8135c375SStefano Zampini                                  0,3,5,4,
311*8135c375SStefano Zampini                                  2,1,7,6,
312*8135c375SStefano Zampini                                  3,2,6,5,
313*8135c375SStefano Zampini                                  0,4,7,1};
314*8135c375SStefano Zampini   PetscContainer       glvis_container;
315*8135c375SStefano Zampini   PetscBool            enabled = PETSC_TRUE;
316*8135c375SStefano Zampini 
317*8135c375SStefano Zampini   PetscFunctionBegin;
318*8135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
319*8135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
320*8135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
321*8135c375SStefano Zampini   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
322*8135c375SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr);
323*8135c375SStefano Zampini   if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
324*8135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
325*8135c375SStefano Zampini   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
326*8135c375SStefano Zampini   if (depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated");
327*8135c375SStefano Zampini 
328*8135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
329*8135c375SStefano Zampini   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
330*8135c375SStefano Zampini   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
331*8135c375SStefano Zampini   {
332*8135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
333*8135c375SStefano Zampini     ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
334*8135c375SStefano Zampini     enabled = glvis_info->enabled;
335*8135c375SStefano Zampini   }
336*8135c375SStefano Zampini   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, &fEndInterior, NULL, NULL);CHKERRQ(ierr);
337*8135c375SStefano Zampini   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
338*8135c375SStefano Zampini   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
339*8135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr);
340*8135c375SStefano Zampini 
341*8135c375SStefano Zampini   /*
342*8135c375SStefano Zampini      a couple of sections of the mesh specification are disabled
343*8135c375SStefano Zampini        - boundary: unless we want to visualize boundary attributes, the boundary is not needed for proper mesh visualization
344*8135c375SStefano Zampini        - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package and be able to derefine the mesh
345*8135c375SStefano Zampini   */
346*8135c375SStefano Zampini   enable_boundary = PETSC_FALSE;
347*8135c375SStefano Zampini   enable_ncmesh   = PETSC_FALSE;
348*8135c375SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr);
349*8135c375SStefano Zampini   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation; useful for debugging purposes",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr);
350*8135c375SStefano Zampini   ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation; useful for debugging purposes",NULL,enable_ncmesh,&enable_ncmesh,NULL);CHKERRQ(ierr);
351*8135c375SStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
352*8135c375SStefano Zampini 
353*8135c375SStefano Zampini   /* Identify possible cells in the overlap */
354*8135c375SStefano Zampini   gNum = NULL;
355*8135c375SStefano Zampini   novl = 0;
356*8135c375SStefano Zampini   ovl  = PETSC_FALSE;
357*8135c375SStefano Zampini   pown = NULL;
358*8135c375SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
359*8135c375SStefano Zampini   if (commsize > 1) {
360*8135c375SStefano Zampini     ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr);
361*8135c375SStefano Zampini     ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr);
362*8135c375SStefano Zampini     for (p=cStart; p<cEnd; p++) {
363*8135c375SStefano Zampini       if (gNum[p-cStart] < 0) {
364*8135c375SStefano Zampini         ovl = PETSC_TRUE;
365*8135c375SStefano Zampini         novl++;
366*8135c375SStefano Zampini       }
367*8135c375SStefano Zampini     }
368*8135c375SStefano Zampini     if (ovl) {
369*8135c375SStefano Zampini       /* it may happen that pown get not destroyed, if the user closes the window while this function is running.
370*8135c375SStefano Zampini          TODO: garbage collector? attach pown to dm?  */
371*8135c375SStefano Zampini       ierr = PetscBTCreate(PetscMax(cEnd-cStart,vEnd-vStart),&pown);CHKERRQ(ierr);
372*8135c375SStefano Zampini     }
373*8135c375SStefano Zampini   }
374*8135c375SStefano Zampini 
375*8135c375SStefano Zampini   if (!enabled) {
376*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
377*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
378*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
379*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
380*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
381*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
382*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
383*8135c375SStefano Zampini     ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
384*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
385*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
386*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
387*8135c375SStefano Zampini     ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
388*8135c375SStefano Zampini     if (globalNum) {
389*8135c375SStefano Zampini       ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
390*8135c375SStefano Zampini     }
391*8135c375SStefano Zampini     PetscFunctionReturn(0);
392*8135c375SStefano Zampini   }
393*8135c375SStefano Zampini 
394*8135c375SStefano Zampini   /* header */
395*8135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
396*8135c375SStefano Zampini 
397*8135c375SStefano Zampini   /* topological dimension */
398*8135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
399*8135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
400*8135c375SStefano Zampini 
401*8135c375SStefano Zampini   /* elements */
402*8135c375SStefano Zampini   /* TODO: is this the label we want to use for marking material IDs?
403*8135c375SStefano Zampini      We should decide to have a single marker for these stuff
404*8135c375SStefano Zampini      Proposal: DMSetMaterialLabel?
405*8135c375SStefano Zampini   */
406*8135c375SStefano Zampini   ierr = DMGetLabel(dm,"Cell sets",&label);CHKERRQ(ierr);
407*8135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
408*8135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr);
409*8135c375SStefano Zampini   for (p=cStart;p<cEnd;p++) {
410*8135c375SStefano Zampini     PetscInt vids[8],i,nv = 0,cid = -1,mid = 1;
411*8135c375SStefano Zampini 
412*8135c375SStefano Zampini     if (ovl) {
413*8135c375SStefano Zampini       if (gNum[p-cStart] < 0) continue;
414*8135c375SStefano Zampini       else {
415*8135c375SStefano Zampini         ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr);
416*8135c375SStefano Zampini       }
417*8135c375SStefano Zampini     }
418*8135c375SStefano Zampini     ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
419*8135c375SStefano Zampini     ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,&nv,vids);CHKERRQ(ierr);
420*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
421*8135c375SStefano Zampini     for (i=0;i<nv;i++) {
422*8135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer," %D",vids[i]);CHKERRQ(ierr);
423*8135c375SStefano Zampini     }
424*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
425*8135c375SStefano Zampini   }
426*8135c375SStefano Zampini 
427*8135c375SStefano Zampini   /* determine orientation of boundary mesh */
428*8135c375SStefano Zampini   if (cEnd-cStart) {
429*8135c375SStefano Zampini     ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr);
430*8135c375SStefano Zampini     switch(dim) {
431*8135c375SStefano Zampini       case 1:
432*8135c375SStefano Zampini         if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc);
433*8135c375SStefano Zampini         faces = fv1;
434*8135c375SStefano Zampini         vpf = 1;
435*8135c375SStefano Zampini         break;
436*8135c375SStefano Zampini       case 2:
437*8135c375SStefano Zampini         switch (fpc) {
438*8135c375SStefano Zampini           case 3:
439*8135c375SStefano Zampini             faces = fv2tri;
440*8135c375SStefano Zampini             vpf   = 2;
441*8135c375SStefano Zampini             break;
442*8135c375SStefano Zampini           case 4:
443*8135c375SStefano Zampini             faces = fv2quad;
444*8135c375SStefano Zampini             vpf   = 2;
445*8135c375SStefano Zampini             break;
446*8135c375SStefano Zampini           default:
447*8135c375SStefano Zampini             SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
448*8135c375SStefano Zampini             break;
449*8135c375SStefano Zampini         }
450*8135c375SStefano Zampini         break;
451*8135c375SStefano Zampini       case 3:
452*8135c375SStefano Zampini         switch (fpc) {
453*8135c375SStefano Zampini           case 4:
454*8135c375SStefano Zampini             faces = fv3tet;
455*8135c375SStefano Zampini             vpf   = 3;
456*8135c375SStefano Zampini             break;
457*8135c375SStefano Zampini           case 6:
458*8135c375SStefano Zampini             faces = fv3hex;
459*8135c375SStefano Zampini             vpf   = 4;
460*8135c375SStefano Zampini             break;
461*8135c375SStefano Zampini           default:
462*8135c375SStefano Zampini             SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc);
463*8135c375SStefano Zampini             break;
464*8135c375SStefano Zampini         }
465*8135c375SStefano Zampini         break;
466*8135c375SStefano Zampini       default:
467*8135c375SStefano Zampini         SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim");
468*8135c375SStefano Zampini         break;
469*8135c375SStefano Zampini     }
470*8135c375SStefano Zampini   }
471*8135c375SStefano Zampini 
472*8135c375SStefano Zampini   /* boundary */
473*8135c375SStefano Zampini   ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr);
474*8135c375SStefano Zampini   fEnd = fEndInterior < 0 ? fEnd : fEndInterior;
475*8135c375SStefano Zampini   if (!ovl) {
476*8135c375SStefano Zampini     for (p=fStart,bf=0;p<fEnd;p++) {
477*8135c375SStefano Zampini       PetscInt supportSize;
478*8135c375SStefano Zampini 
479*8135c375SStefano Zampini       ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
480*8135c375SStefano Zampini       if (supportSize == 1) bf++;
481*8135c375SStefano Zampini     }
482*8135c375SStefano Zampini   } else {
483*8135c375SStefano Zampini     for (p=fStart,bf=0;p<fEnd;p++) {
484*8135c375SStefano Zampini       const PetscInt *support;
485*8135c375SStefano Zampini       PetscInt       i,supportSize;
486*8135c375SStefano Zampini       PetscBool      has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
487*8135c375SStefano Zampini 
488*8135c375SStefano Zampini       ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
489*8135c375SStefano Zampini       ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
490*8135c375SStefano Zampini       for (i=0;i<supportSize;i++) {
491*8135c375SStefano Zampini         if (PetscBTLookup(pown,support[i]-cStart)) has_owned = PETSC_TRUE;
492*8135c375SStefano Zampini         else has_ghost = PETSC_TRUE;
493*8135c375SStefano Zampini       }
494*8135c375SStefano Zampini       if ((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)) bf++;
495*8135c375SStefano Zampini     }
496*8135c375SStefano Zampini   }
497*8135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
498*8135c375SStefano Zampini   if (!enable_boundary) {
499*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
500*8135c375SStefano Zampini   } else {
501*8135c375SStefano Zampini     /* TODO: is this the label we want to use for marking boundary facets?
502*8135c375SStefano Zampini        We should decide to have a single marker for these stuff
503*8135c375SStefano Zampini        Proposal: DMSetBoundaryLabel?
504*8135c375SStefano Zampini     */
505*8135c375SStefano Zampini     ierr = DMGetLabel(dm,"marker",&label);CHKERRQ(ierr);
506*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr);
507*8135c375SStefano Zampini     if (!ovl) {
508*8135c375SStefano Zampini       for (p=fStart;p<fEnd;p++) {
509*8135c375SStefano Zampini         PetscInt supportSize;
510*8135c375SStefano Zampini 
511*8135c375SStefano Zampini         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
512*8135c375SStefano Zampini         if (supportSize == 1) {
513*8135c375SStefano Zampini           const    PetscInt *support, *cone;
514*8135c375SStefano Zampini           PetscInt i,c,v,cid = -1,mid = -1;
515*8135c375SStefano Zampini           PetscInt cellClosureSize,*cellClosure = NULL;
516*8135c375SStefano Zampini 
517*8135c375SStefano Zampini           ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
518*8135c375SStefano Zampini           ierr = DMPlexGetCone(dm,support[0],&cone);CHKERRQ(ierr);
519*8135c375SStefano Zampini           for (c=0;c<fpc;c++)
520*8135c375SStefano Zampini             if (cone[c] == p)
521*8135c375SStefano Zampini               break;
522*8135c375SStefano Zampini 
523*8135c375SStefano Zampini           ierr = DMPlexGetTransitiveClosure(dm,support[0],PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
524*8135c375SStefano Zampini           for (v=0;v<cellClosureSize;v++)
525*8135c375SStefano Zampini             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd)
526*8135c375SStefano Zampini               break;
527*8135c375SStefano Zampini           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
528*8135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
529*8135c375SStefano Zampini           for (i=0;i<vpf;i++) {
530*8135c375SStefano Zampini             ierr = PetscViewerASCIIPrintf(viewer," %D",cellClosure[2*(v+faces[c*vpf+i])]-vStart);CHKERRQ(ierr);
531*8135c375SStefano Zampini           }
532*8135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
533*8135c375SStefano Zampini           ierr = DMPlexRestoreTransitiveClosure(dm,support[0],PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
534*8135c375SStefano Zampini         }
535*8135c375SStefano Zampini       }
536*8135c375SStefano Zampini     } else {
537*8135c375SStefano Zampini       for (p=fStart;p<fEnd;p++) {
538*8135c375SStefano Zampini         const PetscInt *support;
539*8135c375SStefano Zampini         PetscInt       i,supportSize,validcell = -1;
540*8135c375SStefano Zampini         PetscBool      has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE;
541*8135c375SStefano Zampini 
542*8135c375SStefano Zampini         ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr);
543*8135c375SStefano Zampini         ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
544*8135c375SStefano Zampini         for (i=0;i<supportSize;i++) {
545*8135c375SStefano Zampini           if (PetscBTLookup(pown,support[i]-cStart)) {
546*8135c375SStefano Zampini             has_owned = PETSC_TRUE;
547*8135c375SStefano Zampini             validcell = support[i];
548*8135c375SStefano Zampini           } else has_ghost = PETSC_TRUE;
549*8135c375SStefano Zampini         }
550*8135c375SStefano Zampini         if ((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)) {
551*8135c375SStefano Zampini           const    PetscInt *support, *cone;
552*8135c375SStefano Zampini           PetscInt i,c,v,cid = -1,mid = -1;
553*8135c375SStefano Zampini           PetscInt cellClosureSize,*cellClosure = NULL;
554*8135c375SStefano Zampini 
555*8135c375SStefano Zampini           ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr);
556*8135c375SStefano Zampini           ierr = DMPlexGetCone(dm,validcell,&cone);CHKERRQ(ierr);
557*8135c375SStefano Zampini           for (c=0;c<fpc;c++)
558*8135c375SStefano Zampini             if (cone[c] == p)
559*8135c375SStefano Zampini               break;
560*8135c375SStefano Zampini 
561*8135c375SStefano Zampini           ierr = DMPlexGetTransitiveClosure(dm,validcell,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
562*8135c375SStefano Zampini           for (v=0;v<cellClosureSize;v++)
563*8135c375SStefano Zampini             if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd)
564*8135c375SStefano Zampini               break;
565*8135c375SStefano Zampini           ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr);
566*8135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr);
567*8135c375SStefano Zampini           for (i=0;i<vpf;i++) {
568*8135c375SStefano Zampini             ierr = PetscViewerASCIIPrintf(viewer," %D",cellClosure[2*(v+faces[c*vpf+i])]-vStart);CHKERRQ(ierr);
569*8135c375SStefano Zampini           }
570*8135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
571*8135c375SStefano Zampini           ierr = DMPlexRestoreTransitiveClosure(dm,validcell,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr);
572*8135c375SStefano Zampini         }
573*8135c375SStefano Zampini       }
574*8135c375SStefano Zampini     }
575*8135c375SStefano Zampini   }
576*8135c375SStefano Zampini 
577*8135c375SStefano Zampini   /* mark owned vertices */
578*8135c375SStefano Zampini   if (ovl) {
579*8135c375SStefano Zampini     ierr = PetscBTMemzero(vEnd-vStart,pown);CHKERRQ(ierr);
580*8135c375SStefano Zampini     for (p=cStart;p<cEnd;p++) {
581*8135c375SStefano Zampini       PetscInt i,closureSize,*closure = NULL;
582*8135c375SStefano Zampini 
583*8135c375SStefano Zampini       if (gNum[p-cStart] < 0) continue;
584*8135c375SStefano Zampini       ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
585*8135c375SStefano Zampini       for (i=0;i<closureSize;i++) {
586*8135c375SStefano Zampini         const PetscInt pp = closure[2*i];
587*8135c375SStefano Zampini 
588*8135c375SStefano Zampini         if (pp >= vStart && pp < vEnd) {
589*8135c375SStefano Zampini           ierr = PetscBTSet(pown,pp-vStart);CHKERRQ(ierr);
590*8135c375SStefano Zampini         }
591*8135c375SStefano Zampini       }
592*8135c375SStefano Zampini       ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
593*8135c375SStefano Zampini     }
594*8135c375SStefano Zampini   }
595*8135c375SStefano Zampini   if (globalNum) {
596*8135c375SStefano Zampini     ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr);
597*8135c375SStefano Zampini   }
598*8135c375SStefano Zampini 
599*8135c375SStefano Zampini   /* vertex_parents (Non-conforming meshes) */
600*8135c375SStefano Zampini   parentSection  = NULL;
601*8135c375SStefano Zampini   if (enable_ncmesh) {
602*8135c375SStefano Zampini     ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
603*8135c375SStefano Zampini   }
604*8135c375SStefano Zampini   if (parentSection) {
605*8135c375SStefano Zampini     PetscInt vp,gvp;
606*8135c375SStefano Zampini 
607*8135c375SStefano Zampini     for (vp=0,p=vStart;p<vEnd;p++) {
608*8135c375SStefano Zampini       DMLabel  dlabel;
609*8135c375SStefano Zampini       PetscInt parent,depth;
610*8135c375SStefano Zampini 
611*8135c375SStefano Zampini       if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
612*8135c375SStefano Zampini       ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
613*8135c375SStefano Zampini       ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr);
614*8135c375SStefano Zampini       ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
615*8135c375SStefano Zampini       if (parent != p) vp++;
616*8135c375SStefano Zampini     }
617*8135c375SStefano Zampini     ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
618*8135c375SStefano Zampini     if (gvp) {
619*8135c375SStefano Zampini       PetscInt  maxsupp;
620*8135c375SStefano Zampini       PetscBool *skip = NULL;
621*8135c375SStefano Zampini 
622*8135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr);
623*8135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr);
624*8135c375SStefano Zampini       ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr);
625*8135c375SStefano Zampini       ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr);
626*8135c375SStefano Zampini       for (p=vStart;p<vEnd;p++) {
627*8135c375SStefano Zampini         DMLabel  dlabel;
628*8135c375SStefano Zampini         PetscInt parent;
629*8135c375SStefano Zampini 
630*8135c375SStefano Zampini         if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue;
631*8135c375SStefano Zampini         ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr);
632*8135c375SStefano Zampini         ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
633*8135c375SStefano Zampini         if (parent != p) {
634*8135c375SStefano Zampini           PetscInt       vids[8],i,nv,size,n;
635*8135c375SStefano Zampini           PetscInt       numChildren,depth = -1;
636*8135c375SStefano Zampini           const PetscInt *children;
637*8135c375SStefano Zampini           ierr = DMPlexGetConeSize(dm,parent,&size);CHKERRQ(ierr);
638*8135c375SStefano Zampini           switch (size) {
639*8135c375SStefano Zampini             case 2: /* edge */
640*8135c375SStefano Zampini               nv   = 0;
641*8135c375SStefano Zampini               ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,&nv,vids);CHKERRQ(ierr);
642*8135c375SStefano Zampini               ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr);
643*8135c375SStefano Zampini               for (i=0;i<nv;i++) {
644*8135c375SStefano Zampini                 ierr = PetscViewerASCIIPrintf(viewer," %D",vids[i]);CHKERRQ(ierr);
645*8135c375SStefano Zampini               }
646*8135c375SStefano Zampini               ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
647*8135c375SStefano Zampini               vp--;
648*8135c375SStefano Zampini               break;
649*8135c375SStefano Zampini             case 4: /* face */
650*8135c375SStefano Zampini               ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr);
651*8135c375SStefano Zampini               ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
652*8135c375SStefano Zampini               for (n=0;n<numChildren;n++) {
653*8135c375SStefano Zampini                 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr);
654*8135c375SStefano Zampini                 if (!depth) {
655*8135c375SStefano Zampini                   const PetscInt *hvsupp,*hesupp,*cone;
656*8135c375SStefano Zampini                   PetscInt       hvsuppSize,hesuppSize,coneSize;
657*8135c375SStefano Zampini                   PetscInt       hv = children[n],he,f;
658*8135c375SStefano Zampini                   PetscBool      found = PETSC_FALSE;
659*8135c375SStefano Zampini 
660*8135c375SStefano Zampini                   ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr);
661*8135c375SStefano Zampini                   ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr);
662*8135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
663*8135c375SStefano Zampini                     PetscInt ep;
664*8135c375SStefano Zampini                     ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr);
665*8135c375SStefano Zampini                     if (ep != hvsupp[i]) {
666*8135c375SStefano Zampini                       he = hvsupp[i];
667*8135c375SStefano Zampini                       found = PETSC_TRUE;
668*8135c375SStefano Zampini                     } else {
669*8135c375SStefano Zampini                       skip[i] = PETSC_TRUE;
670*8135c375SStefano Zampini                     }
671*8135c375SStefano Zampini                   }
672*8135c375SStefano Zampini                   if (!found) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize);
673*8135c375SStefano Zampini                   he      = hvsupp[i];
674*8135c375SStefano Zampini                   skip[i] = PETSC_TRUE;
675*8135c375SStefano Zampini                   ierr    = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr);
676*8135c375SStefano Zampini                   vids[0] = (cone[0] == hv) ? cone[1] : cone[0];
677*8135c375SStefano Zampini                   ierr    = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr);
678*8135c375SStefano Zampini                   ierr    = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr);
679*8135c375SStefano Zampini                   for (f=0;f<hesuppSize;f++) {
680*8135c375SStefano Zampini                     PetscInt j;
681*8135c375SStefano Zampini 
682*8135c375SStefano Zampini                     ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr);
683*8135c375SStefano Zampini                     ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr);
684*8135c375SStefano Zampini                     for (j=0;j<coneSize;j++) {
685*8135c375SStefano Zampini                       PetscInt k;
686*8135c375SStefano Zampini                       for (k=0;k<hvsuppSize;k++) {
687*8135c375SStefano Zampini                         if (hvsupp[k] == cone[j]) {
688*8135c375SStefano Zampini                           skip[k] = PETSC_TRUE;
689*8135c375SStefano Zampini                           break;
690*8135c375SStefano Zampini                         }
691*8135c375SStefano Zampini                       }
692*8135c375SStefano Zampini                     }
693*8135c375SStefano Zampini                   }
694*8135c375SStefano Zampini                   for (i=0;i<hvsuppSize;i++) {
695*8135c375SStefano Zampini                     if (!skip[i]) {
696*8135c375SStefano Zampini                       ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr);
697*8135c375SStefano Zampini                       vids[1] = (cone[0] == hv) ? cone[1] : cone[0];
698*8135c375SStefano Zampini                     }
699*8135c375SStefano Zampini                   }
700*8135c375SStefano Zampini                   ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr);
701*8135c375SStefano Zampini                   for (i=0;i<2;i++) {
702*8135c375SStefano Zampini                     ierr = PetscViewerASCIIPrintf(viewer," %D",vids[i]-vStart);CHKERRQ(ierr);
703*8135c375SStefano Zampini                   }
704*8135c375SStefano Zampini                   ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
705*8135c375SStefano Zampini                   vp--;
706*8135c375SStefano Zampini                 }
707*8135c375SStefano Zampini               }
708*8135c375SStefano Zampini               break;
709*8135c375SStefano Zampini             default:
710*8135c375SStefano Zampini               SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %d",size);
711*8135c375SStefano Zampini           }
712*8135c375SStefano Zampini         }
713*8135c375SStefano Zampini       }
714*8135c375SStefano Zampini       ierr = PetscFree(skip);CHKERRQ(ierr);
715*8135c375SStefano Zampini     }
716*8135c375SStefano Zampini     if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp);
717*8135c375SStefano Zampini   }
718*8135c375SStefano Zampini   ierr = PetscBTDestroy(&pown);CHKERRQ(ierr);
719*8135c375SStefano Zampini 
720*8135c375SStefano Zampini   /* vertices */
721*8135c375SStefano Zampini   ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
722*8135c375SStefano Zampini   ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr);
723*8135c375SStefano Zampini   ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr);
724*8135c375SStefano Zampini   ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr);
725*8135c375SStefano Zampini   ierr = DMPlexGetDepthStratum(cdm,0,&vStart,&vEnd);CHKERRQ(ierr);
726*8135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
727*8135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr);
728*8135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
729*8135c375SStefano Zampini   ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr);
730*8135c375SStefano Zampini   for (p=vStart;p<vEnd;p++) {
731*8135c375SStefano Zampini     PetscInt i,st,dof;
732*8135c375SStefano Zampini 
733*8135c375SStefano Zampini     ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr);
734*8135c375SStefano Zampini     ierr = PetscSectionGetOffset(coordSection,p,&st);CHKERRQ(ierr);
735*8135c375SStefano Zampini     for (i=st;i<st+dof-1;i++) {
736*8135c375SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"%g ",PetscRealPart(array[i]));CHKERRQ(ierr);
737*8135c375SStefano Zampini     }
738*8135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%g\n",PetscRealPart(array[i]));CHKERRQ(ierr);
739*8135c375SStefano Zampini   }
740*8135c375SStefano Zampini   ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr);
741*8135c375SStefano Zampini   PetscFunctionReturn(0);
742*8135c375SStefano Zampini }
743*8135c375SStefano Zampini 
744*8135c375SStefano Zampini /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */
745*8135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer)
746*8135c375SStefano Zampini {
747*8135c375SStefano Zampini   PetscErrorCode ierr;
748*8135c375SStefano Zampini   PetscBool      isglvis,isascii;
749*8135c375SStefano Zampini 
750*8135c375SStefano Zampini   PetscFunctionBegin;
751*8135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
752*8135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
753*8135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
754*8135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
755*8135c375SStefano Zampini   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
756*8135c375SStefano Zampini   if (isglvis) {
757*8135c375SStefano Zampini     PetscViewer          view;
758*8135c375SStefano Zampini     PetscViewerGLVisType type;
759*8135c375SStefano Zampini 
760*8135c375SStefano Zampini     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
761*8135c375SStefano Zampini     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
762*8135c375SStefano Zampini     if (view) { /* in the socket case, it may happen that the connection failed */
763*8135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
764*8135c375SStefano Zampini         PetscMPIInt size,rank;
765*8135c375SStefano Zampini         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
766*8135c375SStefano Zampini         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
767*8135c375SStefano Zampini         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
768*8135c375SStefano Zampini       }
769*8135c375SStefano Zampini       ierr = DMPlexView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
770*8135c375SStefano Zampini       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
771*8135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
772*8135c375SStefano Zampini         PetscInt    dim;
773*8135c375SStefano Zampini         const char* name;
774*8135c375SStefano Zampini 
775*8135c375SStefano Zampini         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
776*8135c375SStefano Zampini         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
777*8135c375SStefano Zampini         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
778*8135c375SStefano Zampini         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
779*8135c375SStefano Zampini       }
780*8135c375SStefano Zampini     }
781*8135c375SStefano Zampini   } else {
782*8135c375SStefano Zampini     ierr = DMPlexView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
783*8135c375SStefano Zampini   }
784*8135c375SStefano Zampini   PetscFunctionReturn(0);
785*8135c375SStefano Zampini }
786