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 = DMGetLocalSection(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 ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 137 if (continuous && order > 0) { /* no support for high-order viz, still have to figure out the numbering */ 138 ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim);CHKERRQ(ierr); 139 } else { 140 if (!continuous && order) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported for order %D",order); 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 = VecScatterCreate(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 = VecScatterCreate(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 minl, 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 *mid = *mid - minl + 1; /* MFEM does not like negative markers */ 257 } else *mid = 1; 258 if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 259 #if defined PETSC_USE_DEBUG 260 if (dim < 0 || dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim); 261 if (csize > 8) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p); 262 if (depth != 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found depth %D for point %D. You should interpolate the mesh first",depth,p); 263 #endif 264 *cid = mfem_table_cid_unint[dim][csize]; 265 } else { 266 #if defined PETSC_USE_DEBUG 267 if (csize > 6) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p); 268 if (pdepth < 0 || pdepth > 3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p); 269 #endif 270 *cid = mfem_table_cid[pdepth][csize]; 271 } 272 PetscFunctionReturn(0); 273 } 274 275 static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, int vids[]) 276 { 277 PetscInt dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL; 278 PetscErrorCode ierr; 279 280 PetscFunctionBegin; 281 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 282 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 283 sdim = dim; 284 if (csec) { 285 PetscInt sStart,sEnd; 286 287 ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 288 ierr = PetscSectionGetChart(csec,&sStart,&sEnd);CHKERRQ(ierr); 289 ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr); 290 off = off/sdim; 291 if (p >= sStart && p < sEnd) { 292 ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 293 } 294 } 295 if (!dof) { 296 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 297 for (i=0,q=0;i<numPoints*2;i+= 2) 298 if ((points[i] >= vStart) && (points[i] < vEnd)) 299 vids[q++] = (int)(points[i]-vStart+off); 300 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 301 } else { 302 ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr); 303 ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 304 for (q=0;q<dof/sdim;q++) vids[q] = (int)(off/sdim + q); 305 } 306 *nv = q; 307 PetscFunctionReturn(0); 308 } 309 310 /* 311 ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR 312 Higher order meshes are also supported 313 */ 314 static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer) 315 { 316 DMLabel label; 317 PetscSection coordSection,parentSection; 318 Vec coordinates,hovec; 319 const PetscScalar *array; 320 PetscInt bf,p,sdim,dim,depth,novl,minl; 321 PetscInt cStart,cEnd,vStart,vEnd,nvert; 322 PetscMPIInt size; 323 PetscBool localized,isascii; 324 PetscBool enable_mfem,enable_boundary,enable_ncmesh,view_ovl = PETSC_FALSE; 325 PetscBT pown,vown; 326 PetscErrorCode ierr; 327 PetscContainer glvis_container; 328 PetscBool cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE; 329 PetscBool enable_emark,enable_bmark; 330 const char *fmt; 331 char emark[64] = "",bmark[64] = ""; 332 333 PetscFunctionBegin; 334 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 335 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 336 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 337 if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII"); 338 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr); 339 if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); 340 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 341 342 /* get container: determines if a process visualizes is portion of the data or not */ 343 ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr); 344 if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); 345 { 346 PetscViewerGLVisInfo glvis_info; 347 ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr); 348 enabled = glvis_info->enabled; 349 fmt = glvis_info->fmt; 350 } 351 352 /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh 353 DMPlex does not currently support HO meshes, so there's no API for this */ 354 ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr); 355 356 ierr = DMPlexGetSimplexOrBoxCells(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 357 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 358 ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr); 359 ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr); 360 if (periodic && !localized && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized"); 361 ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); 362 ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 363 ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr); 364 if (!coordinates && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector"); 365 366 /* 367 a couple of sections of the mesh specification are disabled 368 - 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) 369 - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package 370 and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel) 371 */ 372 enable_boundary = PETSC_FALSE; 373 enable_ncmesh = PETSC_FALSE; 374 enable_mfem = PETSC_FALSE; 375 enable_emark = PETSC_FALSE; 376 enable_bmark = PETSC_FALSE; 377 /* I'm tired of problems with negative values in the markers, disable them */ 378 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr); 379 ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr); 380 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); 381 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); 382 ierr = PetscOptionsBool("-viewer_glvis_dm_plex_overlap","Include overlap region in local meshes",NULL,view_ovl,&view_ovl,NULL);CHKERRQ(ierr); 383 ierr = PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark);CHKERRQ(ierr); 384 ierr = PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark);CHKERRQ(ierr); 385 ierr = PetscOptionsEnd();CHKERRQ(ierr); 386 if (enable_bmark) enable_boundary = PETSC_TRUE; 387 388 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 389 if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel"); 390 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 391 if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " 392 "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0"); 393 if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " 394 "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0"); 395 if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 396 if (depth != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth); 397 cellvertex = PETSC_TRUE; 398 } 399 400 /* Identify possible cells in the overlap */ 401 novl = 0; 402 pown = NULL; 403 if (size > 1) { 404 IS globalNum = NULL; 405 const PetscInt *gNum; 406 PetscBool ovl = PETSC_FALSE; 407 408 ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr); 409 if (!globalNum) { 410 if (view_ovl) { 411 ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),cEnd-cStart,0,1,&globalNum);CHKERRQ(ierr); 412 } else { 413 ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr); 414 } 415 ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr); 416 ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr); 417 } 418 ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 419 for (p=cStart; p<cEnd; p++) { 420 if (gNum[p-cStart] < 0) { 421 ovl = PETSC_TRUE; 422 novl++; 423 } 424 } 425 if (ovl) { 426 /* it may happen that pown get not destroyed, if the user closes the window while this function is running. 427 TODO: garbage collector? attach pown to dm? */ 428 ierr = PetscBTCreate(cEnd-cStart,&pown);CHKERRQ(ierr); 429 for (p=cStart; p<cEnd; p++) { 430 if (gNum[p-cStart] < 0) continue; 431 else { 432 ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr); 433 } 434 } 435 } 436 ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 437 } 438 439 /* return if this process is disabled */ 440 if (!enabled) { 441 ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 442 ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 443 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 444 ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 445 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 446 ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 447 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 448 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 449 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 450 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 451 ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 452 PetscFunctionReturn(0); 453 } 454 455 if (enable_mfem) { 456 if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */ 457 PetscInt vpc = 0; 458 char fec[64]; 459 int vids[8] = {0,1,2,3,4,5,6,7}; 460 int hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3}; 461 int quadv[8] = {0,1,3,2}, triv[3] = {0,1,2}; 462 int *dof = NULL; 463 PetscScalar *array,*ptr; 464 465 ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim);CHKERRQ(ierr); 466 if (cEnd-cStart) { 467 PetscInt fpc; 468 469 ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 470 switch(dim) { 471 case 1: 472 vpc = 2; 473 dof = hexv; 474 break; 475 case 2: 476 switch (fpc) { 477 case 3: 478 vpc = 3; 479 dof = triv; 480 break; 481 case 4: 482 vpc = 4; 483 dof = quadv; 484 break; 485 default: 486 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 487 break; 488 } 489 break; 490 case 3: 491 switch (fpc) { 492 case 4: /* TODO: still need to understand L2 ordering for tets */ 493 vpc = 4; 494 dof = tetv; 495 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case"); 496 break; 497 case 6: 498 if (cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc); 499 vpc = 8; 500 dof = hexv; 501 break; 502 case 8: 503 if (!cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 504 vpc = 8; 505 dof = hexv; 506 break; 507 default: 508 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 509 break; 510 } 511 break; 512 default: 513 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 514 break; 515 } 516 ierr = DMPlexInvertCell(dim,vpc,vids);CHKERRQ(ierr); 517 } 518 if (!dof) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs"); 519 ierr = VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);CHKERRQ(ierr); 520 ierr = PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)hovec);CHKERRQ(ierr); 521 ierr = PetscObjectDereference((PetscObject)hovec);CHKERRQ(ierr); 522 ierr = PetscObjectSetName((PetscObject)hovec,fec);CHKERRQ(ierr); 523 ierr = VecGetArray(hovec,&array);CHKERRQ(ierr); 524 ptr = array; 525 for (p=cStart;p<cEnd;p++) { 526 PetscInt csize,v,d; 527 PetscScalar *vals = NULL; 528 529 if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 530 ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 531 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); 532 for (v=0;v<vpc;v++) { 533 for (d=0;d<sdim;d++) { 534 ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d]; 535 } 536 } 537 ptr += vpc*sdim; 538 ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 539 } 540 ierr = VecRestoreArray(hovec,&array);CHKERRQ(ierr); 541 } 542 } 543 /* if we have high-order coordinates in 3D, we need to specify the boundary */ 544 if (hovec && dim == 3) enable_boundary = PETSC_TRUE; 545 546 /* header */ 547 ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 548 549 /* topological dimension */ 550 ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 551 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 552 553 /* elements */ 554 minl = 1; 555 label = NULL; 556 if (enable_emark) { 557 PetscInt lminl = PETSC_MAX_INT; 558 559 ierr = DMGetLabel(dm,emark,&label);CHKERRQ(ierr); 560 if (label) { 561 IS vals; 562 PetscInt ldef; 563 564 ierr = DMLabelGetDefaultValue(label,&ldef);CHKERRQ(ierr); 565 ierr = DMLabelGetValueIS(label,&vals);CHKERRQ(ierr); 566 ierr = ISGetMinMax(vals,&lminl,NULL);CHKERRQ(ierr); 567 ierr = ISDestroy(&vals);CHKERRQ(ierr); 568 lminl = PetscMin(ldef,lminl); 569 } 570 ierr = MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 571 if (minl == PETSC_MAX_INT) minl = 1; 572 } 573 ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 574 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr); 575 for (p=cStart;p<cEnd;p++) { 576 int vids[8]; 577 PetscInt i,nv = 0,cid = -1,mid = 1; 578 579 if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 580 ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);CHKERRQ(ierr); 581 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 582 /* TODO Rewrite to invert by cell type */ 583 ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 584 ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 585 for (i=0;i<nv;i++) { 586 ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr); 587 } 588 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 589 } 590 591 /* boundary */ 592 ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 593 if (!enable_boundary) { 594 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 595 } else { 596 DMLabel perLabel; 597 PetscBT bfaces; 598 PetscInt fStart,fEnd,*fcells; 599 PetscInt *faces = NULL,fpc = 0,vpf = 0, vpc = 0; 600 PetscInt *facesH = NULL,fpcH = 0,vpfH = 0; 601 PetscInt fv1[] = {0,1}, 602 fv2tri[] = {0,1, 603 1,2, 604 2,0}, 605 fv2quad[] = {0,1, 606 1,2, 607 2,3, 608 3,0}, 609 fv2quadH[] = {0,1, 610 2,3, 611 0,2, 612 1,3}, 613 fv3tet[] = {0,1,2, 614 0,3,1, 615 0,2,3, 616 2,1,3}, 617 fv3wedge[] = {0,1,2,-1, 618 3,4,5,-1, 619 0,1,3,4, 620 1,2,4,5, 621 2,0,5,3}, 622 fv3hex[] = {0,1,2,3, 623 4,5,6,7, 624 0,3,5,4, 625 2,1,7,6, 626 3,2,6,5, 627 0,4,7,1}; 628 629 /* determine orientation of boundary mesh */ 630 if (cEnd-cStart) { 631 ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 632 switch(dim) { 633 case 1: 634 if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc); 635 faces = fv1; 636 vpf = 1; 637 vpc = 2; 638 break; 639 case 2: 640 switch (fpc) { 641 case 0: 642 case 3: 643 faces = fv2tri; 644 vpf = 2; 645 vpc = 3; 646 if (fpcH == 4) { 647 facesH = fv2quadH; 648 vpfH = 2; 649 } else if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 650 break; 651 case 4: 652 faces = fv2quad; 653 vpf = 2; 654 vpc = 4; 655 if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 656 break; 657 default: 658 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 659 break; 660 } 661 break; 662 case 3: 663 switch (fpc) { 664 case 0: 665 case 4: 666 faces = fv3tet; 667 vpf = 3; 668 vpc = 4; 669 if (fpcH == 5) { 670 facesH = fv3wedge; 671 vpfH = -4; 672 } else if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 673 break; 674 case 6: 675 faces = fv3hex; 676 vpf = 4; 677 vpc = 8; 678 if (fpcH) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case hybrid faces per cell %D",fpcH); 679 break; 680 default: 681 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 682 break; 683 } 684 break; 685 default: 686 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 687 break; 688 } 689 } 690 ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr); 691 ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr); 692 ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr); 693 ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr); 694 ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 695 if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */ 696 ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr); 697 ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 698 ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr); 699 for (p=cStart;p<cEnd;p++) { 700 PetscInt dof, uvpc; 701 702 ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr); 703 if (dof) { 704 PetscInt v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL; 705 PetscScalar *vals = NULL; 706 /* TODO Check cell type */ 707 uvpc = vpc; 708 if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim); 709 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); 710 ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 711 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 712 for (v=0;v<cellClosureSize;v++) 713 if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) { 714 vidxs = cellClosure + 2*v; 715 break; 716 } 717 if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices"); 718 for (v=0;v<uvpc;v++) { 719 PetscInt s; 720 721 for (s=0;s<sdim;s++) { 722 if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+uvpc*sdim])>PETSC_MACHINE_EPSILON) { 723 ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr); 724 } 725 } 726 } 727 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 728 ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 729 } 730 } 731 if (dim > 1) { 732 PetscInt eEnd,eStart; 733 734 ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr); 735 for (p=eStart;p<eEnd;p++) { 736 const PetscInt *cone; 737 PetscInt coneSize,i; 738 PetscBool ispe = PETSC_TRUE; 739 740 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 741 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 742 for (i=0;i<coneSize;i++) { 743 PetscInt v; 744 745 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 746 ispe = (PetscBool)(ispe && (v==2)); 747 } 748 if (ispe && coneSize) { 749 PetscInt ch, numChildren; 750 const PetscInt *children; 751 752 ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 753 ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 754 for (ch = 0; ch < numChildren; ch++) { 755 ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr); 756 } 757 } 758 } 759 if (dim > 2) { 760 for (p=fStart;p<fEnd;p++) { 761 const PetscInt *cone; 762 PetscInt coneSize,i; 763 PetscBool ispe = PETSC_TRUE; 764 765 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 766 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 767 for (i=0;i<coneSize;i++) { 768 PetscInt v; 769 770 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 771 ispe = (PetscBool)(ispe && (v==2)); 772 } 773 if (ispe && coneSize) { 774 PetscInt ch, numChildren; 775 const PetscInt *children; 776 777 ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 778 ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 779 for (ch = 0; ch < numChildren; ch++) { 780 ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr); 781 } 782 } 783 } 784 } 785 } 786 } 787 for (p=fStart;p<fEnd;p++) { 788 const PetscInt *support; 789 PetscInt supportSize; 790 PetscBool isbf = PETSC_FALSE; 791 792 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 793 if (pown) { 794 PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE; 795 PetscInt i; 796 797 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 798 for (i=0;i<supportSize;i++) { 799 if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE; 800 else has_ghost = PETSC_TRUE; 801 } 802 isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)); 803 } else { 804 isbf = (PetscBool)(supportSize == 1); 805 } 806 if (!isbf && perLabel) { 807 const PetscInt *cone; 808 PetscInt coneSize,i; 809 810 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 811 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 812 isbf = PETSC_TRUE; 813 for (i=0;i<coneSize;i++) { 814 PetscInt v,d; 815 816 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 817 ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr); 818 isbf = (PetscBool)(isbf && v != d); 819 } 820 } 821 if (isbf) { 822 ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr); 823 } 824 } 825 /* count boundary faces */ 826 for (p=fStart,bf=0;p<fEnd;p++) { 827 if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 828 const PetscInt *support; 829 PetscInt supportSize,c; 830 831 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 832 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 833 for (c=0;c<supportSize;c++) { 834 const PetscInt *cone; 835 PetscInt cell,cl,coneSize; 836 837 cell = support[c]; 838 if (pown && PetscUnlikely(!PetscBTLookup(pown,cell-cStart))) continue; 839 ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 840 ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 841 for (cl=0;cl<coneSize;cl++) { 842 if (cone[cl] == p) { 843 bf += 1; 844 break; 845 } 846 } 847 } 848 } 849 } 850 minl = 1; 851 label = NULL; 852 if (enable_bmark) { 853 PetscInt lminl = PETSC_MAX_INT; 854 855 ierr = DMGetLabel(dm,bmark,&label);CHKERRQ(ierr); 856 if (label) { 857 IS vals; 858 PetscInt ldef; 859 860 ierr = DMLabelGetDefaultValue(label,&ldef);CHKERRQ(ierr); 861 ierr = DMLabelGetValueIS(label,&vals);CHKERRQ(ierr); 862 ierr = ISGetMinMax(vals,&lminl,NULL);CHKERRQ(ierr); 863 ierr = ISDestroy(&vals);CHKERRQ(ierr); 864 lminl = PetscMin(ldef,lminl); 865 } 866 ierr = MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 867 if (minl == PETSC_MAX_INT) minl = 1; 868 } 869 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr); 870 for (p=fStart;p<fEnd;p++) { 871 if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 872 const PetscInt *support; 873 PetscInt supportSize,c,nc = 0; 874 875 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 876 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 877 if (pown) { 878 for (c=0;c<supportSize;c++) { 879 if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) { 880 fcells[nc++] = support[c]; 881 } 882 } 883 } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c]; 884 for (c=0;c<nc;c++) { 885 const PetscInt *cone; 886 int vids[8]; 887 PetscInt i,coneSize,cell,cl,nv,cid = -1,mid = -1; 888 889 cell = fcells[c]; 890 ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 891 ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 892 for (cl=0;cl<coneSize;cl++) 893 if (cone[cl] == p) 894 break; 895 if (cl == coneSize) continue; 896 897 /* face material id and type */ 898 ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);CHKERRQ(ierr); 899 ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 900 /* vertex ids */ 901 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 902 /* TODO Check celltype of cell */ 903 if (0) { /* Is prism */ 904 PetscInt nv = vpfH, inc = vpfH; 905 if (vpfH < 0) { /* Wedge */ 906 if (cl == 0 || cl == 1) nv = 3; 907 else nv = 4; 908 inc = -vpfH; 909 } 910 for (i=0;i<nv;i++) { 911 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[facesH[cl*inc+i]]);CHKERRQ(ierr); 912 } 913 } else { 914 for (i=0;i<vpf;i++) { 915 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[faces[cl*vpf+i]]);CHKERRQ(ierr); 916 } 917 } 918 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 919 bf -= 1; 920 } 921 } 922 } 923 if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf); 924 ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr); 925 ierr = PetscFree(fcells);CHKERRQ(ierr); 926 } 927 928 /* mark owned vertices */ 929 vown = NULL; 930 if (pown) { 931 ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); 932 for (p=cStart;p<cEnd;p++) { 933 PetscInt i,closureSize,*closure = NULL; 934 935 if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue; 936 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 937 for (i=0;i<closureSize;i++) { 938 const PetscInt pp = closure[2*i]; 939 940 if (pp >= vStart && pp < vEnd) { 941 ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr); 942 } 943 } 944 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 945 } 946 } 947 948 /* vertex_parents (Non-conforming meshes) */ 949 parentSection = NULL; 950 if (enable_ncmesh) { 951 ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 952 } 953 if (parentSection) { 954 PetscInt vp,gvp; 955 956 for (vp=0,p=vStart;p<vEnd;p++) { 957 DMLabel dlabel; 958 PetscInt parent,depth; 959 960 if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 961 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 962 ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr); 963 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 964 if (parent != p) vp++; 965 } 966 ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 967 if (gvp) { 968 PetscInt maxsupp; 969 PetscBool *skip = NULL; 970 971 ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr); 972 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr); 973 ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr); 974 ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr); 975 for (p=vStart;p<vEnd;p++) { 976 DMLabel dlabel; 977 PetscInt parent; 978 979 if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 980 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 981 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 982 if (parent != p) { 983 int vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */ 984 PetscInt i,nv,ssize,n,numChildren,depth = -1; 985 const PetscInt *children; 986 987 ierr = DMPlexGetConeSize(dm,parent,&ssize);CHKERRQ(ierr); 988 switch (ssize) { 989 case 2: /* edge */ 990 nv = 0; 991 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 992 ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 993 ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr); 994 for (i=0;i<nv;i++) { 995 ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)vids[i]);CHKERRQ(ierr); 996 } 997 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 998 vp--; 999 break; 1000 case 4: /* face */ 1001 ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr); 1002 for (n=0;n<numChildren;n++) { 1003 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr); 1004 if (!depth) { 1005 const PetscInt *hvsupp,*hesupp,*cone; 1006 PetscInt hvsuppSize,hesuppSize,coneSize; 1007 PetscInt hv = children[n],he = -1,f; 1008 1009 ierr = PetscArrayzero(skip,maxsupp);CHKERRQ(ierr); 1010 ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr); 1011 ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr); 1012 for (i=0;i<hvsuppSize;i++) { 1013 PetscInt ep; 1014 ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr); 1015 if (ep != hvsupp[i]) { 1016 he = hvsupp[i]; 1017 } else { 1018 skip[i] = PETSC_TRUE; 1019 } 1020 } 1021 if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize); 1022 ierr = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr); 1023 vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 1024 ierr = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr); 1025 ierr = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr); 1026 for (f=0;f<hesuppSize;f++) { 1027 PetscInt j; 1028 1029 ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr); 1030 ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr); 1031 for (j=0;j<coneSize;j++) { 1032 PetscInt k; 1033 for (k=0;k<hvsuppSize;k++) { 1034 if (hvsupp[k] == cone[j]) { 1035 skip[k] = PETSC_TRUE; 1036 break; 1037 } 1038 } 1039 } 1040 } 1041 for (i=0;i<hvsuppSize;i++) { 1042 if (!skip[i]) { 1043 ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr); 1044 vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 1045 } 1046 } 1047 ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr); 1048 for (i=0;i<2;i++) { 1049 ierr = PetscViewerASCIIPrintf(viewer," %D",(PetscInt)(vids[i]-vStart));CHKERRQ(ierr); 1050 } 1051 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1052 vp--; 1053 } 1054 } 1055 break; 1056 default: 1057 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize); 1058 } 1059 } 1060 } 1061 ierr = PetscFree(skip);CHKERRQ(ierr); 1062 } 1063 if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp); 1064 } 1065 ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 1066 ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); 1067 1068 /* vertices */ 1069 if (hovec) { /* higher-order meshes */ 1070 const char *fec; 1071 PetscInt i,n,s; 1072 1073 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 1074 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr); 1075 ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr); 1076 ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr); 1077 ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr); 1078 ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr); 1079 ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr); 1080 ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/ 1081 ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr); 1082 ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr); 1083 if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim); 1084 for (i=0;i<n/sdim;i++) { 1085 for (s=0;s<sdim;s++) { 1086 ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr); 1087 } 1088 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1089 } 1090 ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr); 1091 } else { 1092 ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr); 1093 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 1094 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr); 1095 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 1096 ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr); 1097 for (p=0;p<nvert/sdim;p++) { 1098 PetscInt s; 1099 for (s=0;s<sdim;s++) { 1100 PetscReal v = PetscRealPart(array[p*sdim+s]); 1101 1102 ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscIsInfOrNanReal(v) ? 0.0 : v);CHKERRQ(ierr); 1103 } 1104 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1105 } 1106 ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr); 1107 } 1108 PetscFunctionReturn(0); 1109 } 1110 1111 PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer) 1112 { 1113 PetscErrorCode ierr; 1114 PetscFunctionBegin; 1115 ierr = DMView_GLVis(dm,viewer,DMPlexView_GLVis_ASCII);CHKERRQ(ierr); 1116 PetscFunctionReturn(0); 1117 } 1118