1 #include <petsc/private/glvisviewerimpl.h> 2 #include <petsc/private/petscimpl.h> 3 #include <petsc/private/dmpleximpl.h> 4 #include <petscbt.h> 5 #include <petscdmplex.h> 6 #include <petscsf.h> 7 #include <petscds.h> 8 9 typedef struct { 10 PetscInt nf; 11 VecScatter *scctx; 12 } GLVisViewerCtx; 13 14 static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 15 { 16 GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 17 PetscInt i; 18 PetscErrorCode ierr; 19 20 PetscFunctionBegin; 21 for (i=0;i<ctx->nf;i++) { 22 ierr = VecScatterDestroy(&ctx->scctx[i]);CHKERRQ(ierr); 23 } 24 ierr = PetscFree(ctx->scctx);CHKERRQ(ierr); 25 ierr = PetscFree(vctx);CHKERRQ(ierr); 26 PetscFunctionReturn(0); 27 } 28 29 static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 30 { 31 GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 32 PetscInt f; 33 PetscErrorCode ierr; 34 35 PetscFunctionBegin; 36 for (f=0;f<nf;f++) { 37 ierr = VecScatterBegin(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 38 ierr = VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 39 } 40 PetscFunctionReturn(0); 41 } 42 43 /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */ 44 PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer) 45 { 46 DM dm = (DM)odm; 47 Vec xlocal,xfield,*Ufield; 48 PetscDS ds; 49 IS globalNum,isfield; 50 PetscBT vown; 51 char **fieldname = NULL,**fec_type = NULL; 52 const PetscInt *gNum; 53 PetscInt *nlocal,*bs,*idxs,*dims; 54 PetscInt f,maxfields,nfields,c,totc,totdofs,Nv,cum,i; 55 PetscInt dim,cStart,cEnd,vStart,vEnd; 56 GLVisViewerCtx *ctx; 57 PetscSection s; 58 PetscErrorCode ierr; 59 60 PetscFunctionBegin; 61 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 62 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 63 ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 64 ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr); 65 if (!globalNum) { 66 ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr); 67 ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr); 68 ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr); 69 } 70 ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 71 ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); 72 for (c = cStart, totc = 0; c < cEnd; c++) { 73 if (gNum[c-cStart] >= 0) { 74 PetscInt i,numPoints,*points = NULL; 75 76 totc++; 77 ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 78 for (i=0;i<numPoints*2;i+= 2) { 79 if ((points[i] >= vStart) && (points[i] < vEnd)) { 80 ierr = PetscBTSet(vown,points[i]-vStart);CHKERRQ(ierr); 81 } 82 } 83 ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 84 } 85 } 86 for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscLikely(PetscBTLookup(vown,f))) Nv++; 87 88 ierr = DMCreateLocalVector(dm,&xlocal);CHKERRQ(ierr); 89 ierr = VecGetLocalSize(xlocal,&totdofs);CHKERRQ(ierr); 90 ierr = DMGetSection(dm,&s);CHKERRQ(ierr); 91 ierr = PetscSectionGetNumFields(s,&nfields);CHKERRQ(ierr); 92 for (f=0,maxfields=0;f<nfields;f++) { 93 PetscInt bs; 94 95 ierr = PetscSectionGetFieldComponents(s,f,&bs);CHKERRQ(ierr); 96 maxfields += bs; 97 } 98 ierr = PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield);CHKERRQ(ierr); 99 ierr = PetscNew(&ctx);CHKERRQ(ierr); 100 ierr = PetscCalloc1(maxfields,&ctx->scctx);CHKERRQ(ierr); 101 ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 102 if (ds) { 103 for (f=0;f<nfields;f++) { 104 const char* fname; 105 char name[256]; 106 PetscObject disc; 107 size_t len; 108 109 ierr = PetscSectionGetFieldName(s,f,&fname);CHKERRQ(ierr); 110 ierr = PetscStrlen(fname,&len);CHKERRQ(ierr); 111 if (len) { 112 ierr = PetscStrcpy(name,fname);CHKERRQ(ierr); 113 } else { 114 ierr = PetscSNPrintf(name,256,"Field%D",f);CHKERRQ(ierr); 115 } 116 ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr); 117 if (disc) { 118 PetscClassId id; 119 PetscInt Nc; 120 char fec[64]; 121 122 ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr); 123 if (id == PETSCFE_CLASSID) { 124 PetscFE fem = (PetscFE)disc; 125 PetscDualSpace sp; 126 PetscDualSpaceType spname; 127 PetscInt order; 128 PetscBool islag,continuous,H1 = PETSC_TRUE; 129 130 ierr = PetscFEGetNumComponents(fem,&Nc);CHKERRQ(ierr); 131 ierr = PetscFEGetDualSpace(fem,&sp);CHKERRQ(ierr); 132 ierr = PetscDualSpaceGetType(sp,&spname);CHKERRQ(ierr); 133 ierr = PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);CHKERRQ(ierr); 134 if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space"); 135 ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr); 136 if (!continuous) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported"); 137 ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 138 if (continuous && order > 0) { 139 ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim);CHKERRQ(ierr); 140 } else { 141 H1 = PETSC_FALSE; 142 ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P%D",dim,order);CHKERRQ(ierr); 143 } 144 ierr = PetscStrallocpy(name,&fieldname[ctx->nf]);CHKERRQ(ierr); 145 bs[ctx->nf] = Nc; 146 dims[ctx->nf] = dim; 147 if (H1) { 148 nlocal[ctx->nf] = Nc * Nv; 149 ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 150 ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr); 151 for (i=0,cum=0;i<vEnd-vStart;i++) { 152 PetscInt j,off; 153 154 if (PetscUnlikely(!PetscBTLookup(vown,i))) continue; 155 ierr = PetscSectionGetFieldOffset(s,i+vStart,f,&off);CHKERRQ(ierr); 156 for (j=0;j<Nc;j++) idxs[cum++] = off + j; 157 } 158 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 159 } else { 160 nlocal[ctx->nf] = Nc * totc; 161 ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 162 ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr); 163 for (i=0,cum=0;i<cEnd-cStart;i++) { 164 PetscInt j,off; 165 166 if (PetscUnlikely(gNum[i] < 0)) continue; 167 ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr); 168 for (j=0;j<Nc;j++) idxs[cum++] = off + j; 169 } 170 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 171 } 172 ierr = VecScatterCreateWithData(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr); 173 ierr = VecDestroy(&xfield);CHKERRQ(ierr); 174 ierr = ISDestroy(&isfield);CHKERRQ(ierr); 175 ctx->nf++; 176 } else if (id == PETSCFV_CLASSID) { 177 PetscInt c; 178 179 ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr); 180 ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P0",dim);CHKERRQ(ierr); 181 for (c = 0; c < Nc; c++) { 182 char comp[256]; 183 ierr = PetscSNPrintf(comp,256,"%s-Comp%D",name,c);CHKERRQ(ierr); 184 ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr); 185 bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */ 186 nlocal[ctx->nf] = totc; 187 dims[ctx->nf] = dim; 188 ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 189 ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr); 190 for (i=0,cum=0;i<cEnd-cStart;i++) { 191 PetscInt off; 192 193 if (PetscUnlikely(gNum[i])<0) continue; 194 ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr); 195 idxs[cum++] = off + c; 196 } 197 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 198 ierr = VecScatterCreateWithData(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr); 199 ierr = VecDestroy(&xfield);CHKERRQ(ierr); 200 ierr = ISDestroy(&isfield);CHKERRQ(ierr); 201 ctx->nf++; 202 } 203 } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f); 204 } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f); 205 } 206 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM"); 207 ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); 208 ierr = VecDestroy(&xlocal);CHKERRQ(ierr); 209 ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 210 211 /* create work vectors */ 212 for (f=0;f<ctx->nf;f++) { 213 ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]);CHKERRQ(ierr); 214 ierr = PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]);CHKERRQ(ierr); 215 ierr = VecSetBlockSize(Ufield[f],bs[f]);CHKERRQ(ierr); 216 ierr = VecSetDM(Ufield[f],dm);CHKERRQ(ierr); 217 } 218 219 /* customize the viewer */ 220 ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr); 221 for (f=0;f<ctx->nf;f++) { 222 ierr = PetscFree(fieldname[f]);CHKERRQ(ierr); 223 ierr = PetscFree(fec_type[f]);CHKERRQ(ierr); 224 ierr = VecDestroy(&Ufield[f]);CHKERRQ(ierr); 225 } 226 ierr = PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield);CHKERRQ(ierr); 227 PetscFunctionReturn(0); 228 } 229 230 typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_PRISM,MFEM_UNDEF} MFEM_cid; 231 232 MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_UNDEF}, 233 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_UNDEF}, 234 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF,MFEM_UNDEF}, 235 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_PRISM,MFEM_CUBE } }; 236 237 MFEM_cid mfem_table_cid_unint[4][9] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, 238 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, 239 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, 240 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_CUBE } }; 241 242 static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt p, PetscInt *mid, PetscInt *cid) 243 { 244 DMLabel dlabel; 245 PetscInt depth,csize,pdepth,dim; 246 PetscErrorCode ierr; 247 248 PetscFunctionBegin; 249 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 250 ierr = DMLabelGetValue(dlabel,p,&pdepth);CHKERRQ(ierr); 251 ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr); 252 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 253 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 254 if (label) { 255 ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr); 256 } else *mid = 1; 257 if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 258 #if defined PETSC_USE_DEBUG 259 if (dim < 0 || dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim); 260 if (csize > 8) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p); 261 if (depth != 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found depth %D for point %D. You should interpolate the mesh first",depth,p); 262 #endif 263 *cid = mfem_table_cid_unint[dim][csize]; 264 } else { 265 #if defined PETSC_USE_DEBUG 266 if (csize > 6) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p); 267 if (pdepth < 0 || pdepth > 3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p); 268 #endif 269 *cid = mfem_table_cid[pdepth][csize]; 270 } 271 PetscFunctionReturn(0); 272 } 273 274 static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, int vids[]) 275 { 276 PetscInt dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL; 277 PetscErrorCode ierr; 278 279 PetscFunctionBegin; 280 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 281 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 282 sdim = dim; 283 if (csec) { 284 PetscInt sStart,sEnd; 285 286 ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 287 ierr = PetscSectionGetChart(csec,&sStart,&sEnd);CHKERRQ(ierr); 288 ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr); 289 off = off/sdim; 290 if (p >= sStart && p < sEnd) { 291 ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 292 } 293 } 294 if (!dof) { 295 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 296 for (i=0,q=0;i<numPoints*2;i+= 2) 297 if ((points[i] >= vStart) && (points[i] < vEnd)) 298 vids[q++] = (int)(points[i]-vStart+off); 299 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 300 } else { 301 ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr); 302 ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 303 for (q=0;q<dof/sdim;q++) vids[q] = (int)(off/sdim + q); 304 } 305 *nv = q; 306 PetscFunctionReturn(0); 307 } 308 309 static PetscErrorCode DMPlexGlvisInvertHybrid(PetscInt cid,int vids[]) 310 { 311 int tmp; 312 313 PetscFunctionBegin; 314 if (cid == MFEM_SQUARE) { /* PETSc stores hybrid quads not as counter-clockwise quad */ 315 tmp = vids[2]; 316 vids[2] = vids[3]; 317 vids[3] = tmp; 318 } else if (cid == MFEM_PRISM) { /* MFEM uses a different orientation for the base and top triangles of the wedge */ 319 tmp = vids[1]; 320 vids[1] = vids[2]; 321 vids[2] = tmp; 322 tmp = vids[4]; 323 vids[4] = vids[5]; 324 vids[5] = tmp; 325 } 326 PetscFunctionReturn(0); 327 } 328 329 /* 330 ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR 331 Higher order meshes are also supported 332 */ 333 static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer) 334 { 335 DMLabel label; 336 PetscSection coordSection,parentSection; 337 Vec coordinates,hovec; 338 const PetscScalar *array; 339 PetscInt bf,p,sdim,dim,depth,novl; 340 PetscInt cStart,cEnd,cEndInterior,vStart,vEnd,nvert; 341 PetscMPIInt size; 342 PetscBool localized,isascii; 343 PetscBool enable_mfem,enable_boundary,enable_ncmesh; 344 PetscBT pown,vown; 345 PetscErrorCode ierr; 346 PetscContainer glvis_container; 347 PetscBool cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE; 348 const char *fmt; 349 char emark[64] = "",bmark[64] = ""; 350 351 PetscFunctionBegin; 352 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 353 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 354 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 355 if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII"); 356 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr); 357 if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); 358 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 359 360 /* get container: determines if a process visualizes is portion of the data or not */ 361 ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr); 362 if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); 363 { 364 PetscViewerGLVisInfo glvis_info; 365 ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr); 366 enabled = glvis_info->enabled; 367 fmt = glvis_info->fmt; 368 } 369 370 /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh 371 DMPlex does not currently support HO meshes, so there's no API for this */ 372 ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr); 373 374 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 375 ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 376 cEndInterior = cEndInterior < 0 ? cEnd : cEndInterior; 377 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 378 ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr); 379 ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr); 380 if (periodic && !localized && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized"); 381 ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); 382 ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 383 ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr); 384 if (!coordinates && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector"); 385 386 /* 387 a couple of sections of the mesh specification are disabled 388 - boundary: the boundary is not needed for proper mesh visualization unless we want to visualize boundary attributes or we have high-order coordinates in 3D (topologically) 389 - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package 390 and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel) 391 */ 392 enable_boundary = PETSC_FALSE; 393 enable_ncmesh = PETSC_FALSE; 394 enable_mfem = PETSC_FALSE; 395 /* I'm tired of problems with negative values in the markers, disable them */ 396 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr); 397 ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr); 398 ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation (allows derefinement)",NULL,enable_ncmesh,&enable_ncmesh,NULL);CHKERRQ(ierr); 399 ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_mfem","Dump a mesh that can be used with MFEM's FiniteElementSpaces",NULL,enable_mfem,&enable_mfem,NULL);CHKERRQ(ierr); 400 ierr = PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),NULL);CHKERRQ(ierr); 401 ierr = PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),NULL);CHKERRQ(ierr); 402 ierr = PetscOptionsEnd();CHKERRQ(ierr); 403 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 404 if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel"); 405 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 406 if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " 407 "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0"); 408 if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " 409 "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0"); 410 if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 411 if (depth != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth); 412 cellvertex = PETSC_TRUE; 413 } 414 415 /* Identify possible cells in the overlap */ 416 novl = 0; 417 pown = NULL; 418 if (size > 1) { 419 IS globalNum = NULL; 420 const PetscInt *gNum; 421 PetscBool ovl = PETSC_FALSE; 422 423 ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr); 424 if (!globalNum) { 425 ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr); 426 ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr); 427 ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr); 428 } 429 ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 430 for (p=cStart; p<cEnd; p++) { 431 if (gNum[p-cStart] < 0) { 432 ovl = PETSC_TRUE; 433 novl++; 434 } 435 } 436 if (ovl) { 437 /* it may happen that pown get not destroyed, if the user closes the window while this function is running. 438 TODO: garbage collector? attach pown to dm? */ 439 ierr = PetscBTCreate(cEnd-cStart,&pown);CHKERRQ(ierr); 440 for (p=cStart; p<cEnd; p++) { 441 if (gNum[p-cStart] < 0) continue; 442 else { 443 ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr); 444 } 445 } 446 } 447 ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 448 } 449 450 /* return if this process is disabled */ 451 if (!enabled) { 452 ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 453 ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 454 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 455 ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 456 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 457 ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 458 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 459 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 460 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 461 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 462 ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 463 PetscFunctionReturn(0); 464 } 465 466 if (enable_mfem) { 467 if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */ 468 PetscInt vpc = 0; 469 char fec[64]; 470 int vids[8] = {0,1,2,3,4,5,6,7}; 471 int hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3}; 472 int quadv[8] = {0,1,3,2}, triv[3] = {0,1,2}; 473 int *dof = NULL; 474 PetscScalar *array,*ptr; 475 476 ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim);CHKERRQ(ierr); 477 if (cEndInterior < cEnd) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Support for hybrid meshed not currently implemented"); 478 if (cEnd-cStart) { 479 PetscInt fpc; 480 481 ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 482 switch(dim) { 483 case 1: 484 vpc = 2; 485 dof = hexv; 486 break; 487 case 2: 488 switch (fpc) { 489 case 3: 490 vpc = 3; 491 dof = triv; 492 break; 493 case 4: 494 vpc = 4; 495 dof = quadv; 496 break; 497 default: 498 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 499 break; 500 } 501 break; 502 case 3: 503 switch (fpc) { 504 case 4: /* TODO: still need to understand L2 ordering for tets */ 505 vpc = 4; 506 dof = tetv; 507 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case"); 508 break; 509 case 6: 510 if (cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc); 511 vpc = 8; 512 dof = hexv; 513 break; 514 case 8: 515 if (!cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 516 vpc = 8; 517 dof = hexv; 518 break; 519 default: 520 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 521 break; 522 } 523 break; 524 default: 525 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 526 break; 527 } 528 ierr = DMPlexInvertCell(dim,vpc,vids);CHKERRQ(ierr); 529 } 530 if (!dof) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs"); 531 ierr = VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);CHKERRQ(ierr); 532 ierr = PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)hovec);CHKERRQ(ierr); 533 ierr = PetscObjectDereference((PetscObject)hovec);CHKERRQ(ierr); 534 ierr = PetscObjectSetName((PetscObject)hovec,fec);CHKERRQ(ierr); 535 ierr = VecGetArray(hovec,&array);CHKERRQ(ierr); 536 ptr = array; 537 for (p=cStart;p<cEnd;p++) { 538 PetscInt csize,v,d; 539 PetscScalar *vals = NULL; 540 541 if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 542 ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 543 if (csize != vpc*sdim && csize != vpc*sdim*2) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported closure size %D (vpc %D, sdim %D)",csize,vpc,sdim); 544 for (v=0;v<vpc;v++) { 545 for (d=0;d<sdim;d++) { 546 ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d]; 547 } 548 } 549 ptr += vpc*sdim; 550 ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 551 } 552 ierr = VecRestoreArray(hovec,&array);CHKERRQ(ierr); 553 } 554 } 555 /* if we have high-order coordinates in 3D, we need to specify the boundary */ 556 if (hovec && dim == 3) enable_boundary = PETSC_TRUE; 557 558 /* header */ 559 ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 560 561 /* topological dimension */ 562 ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 563 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 564 565 /* elements */ 566 ierr = DMGetLabel(dm,emark,&label);CHKERRQ(ierr); 567 ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 568 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr); 569 for (p=cStart;p<cEnd;p++) { 570 int vids[8]; 571 PetscInt i,nv = 0,cid = -1,mid = 1; 572 573 if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 574 ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr); 575 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 576 ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 577 if (p >= cEndInterior) { 578 ierr = DMPlexGlvisInvertHybrid(cid,vids);CHKERRQ(ierr); 579 } 580 ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 581 for (i=0;i<nv;i++) { 582 ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr); 583 } 584 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 585 } 586 587 /* boundary */ 588 ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 589 if (!enable_boundary) { 590 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 591 } else { 592 DMLabel perLabel; 593 PetscBT bfaces; 594 PetscInt fStart,fEnd,*fcells; 595 PetscInt *faces = NULL,fpc = 0,vpf = 0, vpc = 0; 596 PetscInt *facesH = NULL,fpcH = 0,vpfH = 0, vpcH = 0; 597 PetscInt fv1[] = {0,1}, 598 fv2tri[] = {0,1, 599 1,2, 600 2,0}, 601 fv2quad[] = {0,1, 602 1,2, 603 2,3, 604 3,0}, 605 fv2quadH[] = {0,1, 606 2,3, 607 0,2, 608 1,3}, 609 fv3tet[] = {0,1,2, 610 0,3,1, 611 0,2,3, 612 2,1,3}, 613 fv3wedge[] = {0,1,2,-1, 614 3,4,5,-1, 615 0,1,3,4, 616 1,2,4,5, 617 2,0,5,3}, 618 fv3hex[] = {0,1,2,3, 619 4,5,6,7, 620 0,3,5,4, 621 2,1,7,6, 622 3,2,6,5, 623 0,4,7,1}; 624 625 /* determine orientation of boundary mesh */ 626 if (cEnd-cStart) { 627 if (cEndInterior < cEnd) { 628 ierr = DMPlexGetConeSize(dm,cEndInterior,&fpcH);CHKERRQ(ierr); 629 } 630 if (cEndInterior > cStart) { 631 ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 632 } 633 switch(dim) { 634 case 1: 635 if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc); 636 faces = fv1; 637 vpf = 1; 638 vpc = 2; 639 break; 640 case 2: 641 switch (fpc) { 642 case 0: 643 case 3: 644 faces = fv2tri; 645 vpf = 2; 646 vpc = 3; 647 if (fpcH == 4) { 648 facesH = fv2quadH; 649 vpfH = 2; 650 vpcH = 4; 651 } else if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 652 break; 653 case 4: 654 faces = fv2quad; 655 vpf = 2; 656 vpc = 4; 657 if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 658 break; 659 default: 660 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 661 break; 662 } 663 break; 664 case 3: 665 switch (fpc) { 666 case 0: 667 case 4: 668 faces = fv3tet; 669 vpf = 3; 670 vpc = 4; 671 if (fpcH == 5) { 672 facesH = fv3wedge; 673 vpfH = -4; 674 vpcH = 6; 675 } else if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 676 break; 677 case 6: 678 faces = fv3hex; 679 vpf = 4; 680 vpc = 8; 681 if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 682 break; 683 default: 684 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 685 break; 686 } 687 break; 688 default: 689 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 690 break; 691 } 692 } 693 ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr); 694 ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr); 695 ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr); 696 ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr); 697 ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 698 if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */ 699 ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr); 700 ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 701 ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr); 702 for (p=cStart;p<cEnd;p++) { 703 PetscInt dof, uvpc; 704 705 ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr); 706 if (dof) { 707 PetscInt v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL; 708 PetscScalar *vals = NULL; 709 if (p < cEndInterior) uvpc = vpc; 710 else uvpc = vpcH; 711 if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim); 712 if (dof/sdim != uvpc) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible number of cell dofs %D, vertices %D and space dim %D",dof/sdim,uvpc,sdim); 713 ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 714 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 715 for (v=0;v<cellClosureSize;v++) 716 if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) { 717 vidxs = cellClosure + 2*v; 718 break; 719 } 720 if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices"); 721 for (v=0;v<uvpc;v++) { 722 PetscInt s; 723 724 for (s=0;s<sdim;s++) { 725 if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+uvpc*sdim])>PETSC_MACHINE_EPSILON) { 726 ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr); 727 } 728 } 729 } 730 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 731 ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 732 } 733 } 734 if (dim > 1) { 735 PetscInt eEnd,eStart; 736 737 ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr); 738 for (p=eStart;p<eEnd;p++) { 739 const PetscInt *cone; 740 PetscInt coneSize,i; 741 PetscBool ispe = PETSC_TRUE; 742 743 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 744 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 745 for (i=0;i<coneSize;i++) { 746 PetscInt v; 747 748 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 749 ispe = (PetscBool)(ispe && (v==2)); 750 } 751 if (ispe && coneSize) { 752 PetscInt ch, numChildren; 753 const PetscInt *children; 754 755 ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 756 ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 757 for (ch = 0; ch < numChildren; ch++) { 758 ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr); 759 } 760 } 761 } 762 if (dim > 2) { 763 for (p=fStart;p<fEnd;p++) { 764 const PetscInt *cone; 765 PetscInt coneSize,i; 766 PetscBool ispe = PETSC_TRUE; 767 768 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 769 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 770 for (i=0;i<coneSize;i++) { 771 PetscInt v; 772 773 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 774 ispe = (PetscBool)(ispe && (v==2)); 775 } 776 if (ispe && coneSize) { 777 PetscInt ch, numChildren; 778 const PetscInt *children; 779 780 ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 781 ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 782 for (ch = 0; ch < numChildren; ch++) { 783 ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr); 784 } 785 } 786 } 787 } 788 } 789 } 790 for (p=fStart;p<fEnd;p++) { 791 const PetscInt *support; 792 PetscInt supportSize; 793 PetscBool isbf = PETSC_FALSE; 794 795 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 796 if (pown) { 797 PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE; 798 PetscInt i; 799 800 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 801 for (i=0;i<supportSize;i++) { 802 if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE; 803 else has_ghost = PETSC_TRUE; 804 } 805 isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)); 806 } else { 807 isbf = (PetscBool)(supportSize == 1); 808 } 809 if (!isbf && perLabel) { 810 const PetscInt *cone; 811 PetscInt coneSize,i; 812 813 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 814 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 815 isbf = PETSC_TRUE; 816 for (i=0;i<coneSize;i++) { 817 PetscInt v,d; 818 819 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 820 ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr); 821 isbf = (PetscBool)(isbf && v != d); 822 } 823 } 824 if (isbf) { 825 ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr); 826 } 827 } 828 /* count boundary faces */ 829 for (p=fStart,bf=0;p<fEnd;p++) { 830 if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 831 const PetscInt *support; 832 PetscInt supportSize,c; 833 834 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 835 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 836 for (c=0;c<supportSize;c++) { 837 const PetscInt *cone; 838 PetscInt cell,cl,coneSize; 839 840 cell = support[c]; 841 if (pown && PetscUnlikely(!PetscBTLookup(pown,cell-cStart))) continue; 842 ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 843 ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 844 for (cl=0;cl<coneSize;cl++) { 845 if (cone[cl] == p) { 846 bf += 1; 847 break; 848 } 849 } 850 } 851 } 852 } 853 ierr = DMGetLabel(dm,bmark,&label);CHKERRQ(ierr); 854 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr); 855 for (p=fStart;p<fEnd;p++) { 856 if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 857 const PetscInt *support; 858 PetscInt supportSize,c,nc = 0; 859 860 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 861 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 862 if (pown) { 863 for (c=0;c<supportSize;c++) { 864 if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) { 865 fcells[nc++] = support[c]; 866 } 867 } 868 } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c]; 869 for (c=0;c<nc;c++) { 870 const PetscInt *cone; 871 int vids[8]; 872 PetscInt i,coneSize,cell,cl,nv,cid = -1,mid = -1; 873 874 cell = fcells[c]; 875 ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 876 ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 877 for (cl=0;cl<coneSize;cl++) 878 if (cone[cl] == p) 879 break; 880 if (cl == coneSize) continue; 881 882 /* face material id and type */ 883 ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr); 884 ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 885 /* vertex ids */ 886 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 887 if (cell >= cEndInterior) { 888 PetscInt nv = vpfH, inc = vpfH; 889 if (vpfH < 0) { /* Wedge */ 890 if (cl == 0 || cl == 1) nv = 3; 891 else nv = 4; 892 inc = -vpfH; 893 } 894 for (i=0;i<nv;i++) { 895 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[facesH[cl*inc+i]]);CHKERRQ(ierr); 896 } 897 } else { 898 for (i=0;i<vpf;i++) { 899 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[faces[cl*vpf+i]]);CHKERRQ(ierr); 900 } 901 } 902 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 903 bf -= 1; 904 } 905 } 906 } 907 if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf); 908 ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr); 909 ierr = PetscFree(fcells);CHKERRQ(ierr); 910 } 911 912 /* mark owned vertices */ 913 vown = NULL; 914 if (pown) { 915 ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); 916 for (p=cStart;p<cEnd;p++) { 917 PetscInt i,closureSize,*closure = NULL; 918 919 if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue; 920 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 921 for (i=0;i<closureSize;i++) { 922 const PetscInt pp = closure[2*i]; 923 924 if (pp >= vStart && pp < vEnd) { 925 ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr); 926 } 927 } 928 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 929 } 930 } 931 932 /* vertex_parents (Non-conforming meshes) */ 933 parentSection = NULL; 934 if (enable_ncmesh) { 935 ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 936 } 937 if (parentSection) { 938 PetscInt vp,gvp; 939 940 for (vp=0,p=vStart;p<vEnd;p++) { 941 DMLabel dlabel; 942 PetscInt parent,depth; 943 944 if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 945 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 946 ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr); 947 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 948 if (parent != p) vp++; 949 } 950 ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 951 if (gvp) { 952 PetscInt maxsupp; 953 PetscBool *skip = NULL; 954 955 ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr); 956 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr); 957 ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr); 958 ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr); 959 for (p=vStart;p<vEnd;p++) { 960 DMLabel dlabel; 961 PetscInt parent; 962 963 if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 964 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 965 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 966 if (parent != p) { 967 int vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */ 968 PetscInt i,nv,ssize,n,numChildren,depth = -1; 969 const PetscInt *children; 970 971 ierr = DMPlexGetConeSize(dm,parent,&ssize);CHKERRQ(ierr); 972 switch (ssize) { 973 case 2: /* edge */ 974 nv = 0; 975 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 976 ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 977 ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr); 978 for (i=0;i<nv;i++) { 979 ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr); 980 } 981 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 982 vp--; 983 break; 984 case 4: /* face */ 985 ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr); 986 for (n=0;n<numChildren;n++) { 987 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr); 988 if (!depth) { 989 const PetscInt *hvsupp,*hesupp,*cone; 990 PetscInt hvsuppSize,hesuppSize,coneSize; 991 PetscInt hv = children[n],he = -1,f; 992 993 ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr); 994 ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr); 995 ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr); 996 for (i=0;i<hvsuppSize;i++) { 997 PetscInt ep; 998 ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr); 999 if (ep != hvsupp[i]) { 1000 he = hvsupp[i]; 1001 } else { 1002 skip[i] = PETSC_TRUE; 1003 } 1004 } 1005 if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize); 1006 ierr = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr); 1007 vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 1008 ierr = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr); 1009 ierr = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr); 1010 for (f=0;f<hesuppSize;f++) { 1011 PetscInt j; 1012 1013 ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr); 1014 ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr); 1015 for (j=0;j<coneSize;j++) { 1016 PetscInt k; 1017 for (k=0;k<hvsuppSize;k++) { 1018 if (hvsupp[k] == cone[j]) { 1019 skip[k] = PETSC_TRUE; 1020 break; 1021 } 1022 } 1023 } 1024 } 1025 for (i=0;i<hvsuppSize;i++) { 1026 if (!skip[i]) { 1027 ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr); 1028 vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 1029 } 1030 } 1031 ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr); 1032 for (i=0;i<2;i++) { 1033 ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)(vids[i]-vStart));CHKERRQ(ierr); 1034 } 1035 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1036 vp--; 1037 } 1038 } 1039 break; 1040 default: 1041 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize); 1042 } 1043 } 1044 } 1045 ierr = PetscFree(skip);CHKERRQ(ierr); 1046 } 1047 if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp); 1048 } 1049 ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 1050 ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); 1051 1052 /* vertices */ 1053 if (hovec) { /* higher-order meshes */ 1054 const char *fec; 1055 PetscInt i,n,s; 1056 1057 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 1058 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr); 1059 ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr); 1060 ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr); 1061 ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr); 1062 ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr); 1063 ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr); 1064 ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/ 1065 ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr); 1066 ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr); 1067 if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim); 1068 for (i=0;i<n/sdim;i++) { 1069 for (s=0;s<sdim;s++) { 1070 ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr); 1071 } 1072 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1073 } 1074 ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr); 1075 } else { 1076 ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr); 1077 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 1078 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr); 1079 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 1080 ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr); 1081 for (p=0;p<nvert/sdim;p++) { 1082 PetscInt s; 1083 for (s=0;s<sdim;s++) { 1084 PetscReal v = PetscRealPart(array[p*sdim+s]); 1085 1086 ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscIsInfOrNanReal(v) ? 0.0 : v);CHKERRQ(ierr); 1087 } 1088 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1089 } 1090 ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr); 1091 } 1092 PetscFunctionReturn(0); 1093 } 1094 1095 PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer) 1096 { 1097 PetscErrorCode ierr; 1098 PetscFunctionBegin; 1099 ierr = DMView_GLVis(dm,viewer,DMPlexView_GLVis_ASCII);CHKERRQ(ierr); 1100 PetscFunctionReturn(0); 1101 } 1102