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