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,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 (PetscLikely(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 = PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield);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 /* create work vectors */ 209 for (f=0;f<ctx->nf;f++) { 210 ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]);CHKERRQ(ierr); 211 ierr = PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]);CHKERRQ(ierr); 212 ierr = VecSetBlockSize(Ufield[f],bs[f]);CHKERRQ(ierr); 213 ierr = VecSetDM(Ufield[f],dm);CHKERRQ(ierr); 214 } 215 216 /* customize the viewer */ 217 ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr); 218 for (f=0;f<ctx->nf;f++) { 219 ierr = PetscFree(fieldname[f]);CHKERRQ(ierr); 220 ierr = PetscFree(fec_type[f]);CHKERRQ(ierr); 221 ierr = VecDestroy(&Ufield[f]);CHKERRQ(ierr); 222 } 223 ierr = PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield);CHKERRQ(ierr); 224 PetscFunctionReturn(0); 225 } 226 227 typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_UNDEF} MFEM_cid; 228 229 MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF, MFEM_UNDEF}, 230 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF, MFEM_UNDEF}, 231 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF, MFEM_UNDEF}, 232 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_UNDEF, MFEM_CUBE } }; 233 234 static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt p, PetscInt *mid, PetscInt *cid) 235 { 236 DMLabel dlabel; 237 PetscInt depth,csize; 238 PetscErrorCode ierr; 239 240 PetscFunctionBegin; 241 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 242 ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr); 243 ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr); 244 if (label) { 245 ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr); 246 } else *mid = 1; 247 *cid = mfem_table_cid[depth][csize]; 248 PetscFunctionReturn(0); 249 } 250 251 static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, int vids[]) 252 { 253 PetscInt dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL; 254 PetscErrorCode ierr; 255 256 PetscFunctionBegin; 257 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 258 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 259 sdim = dim; 260 if (csec) { 261 ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 262 ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr); 263 off = off/sdim; 264 ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 265 } 266 if (!dof) { 267 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 268 for (i=0,q=0;i<numPoints*2;i+= 2) 269 if ((points[i] >= vStart) && (points[i] < vEnd)) 270 vids[q++] = points[i]-vStart+off; 271 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 272 } else { 273 ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr); 274 ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 275 for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q; 276 } 277 *nv = q; 278 PetscFunctionReturn(0); 279 } 280 281 /* 282 ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR 283 Higher order meshes are also supported 284 */ 285 static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer) 286 { 287 DMLabel label; 288 PetscSection coordSection,parentSection; 289 Vec coordinates,hovec; 290 const PetscScalar *array; 291 PetscInt bf,p,sdim,dim,depth,novl; 292 PetscInt cStart,cEnd,cEndInterior,vStart,vEnd,nvert; 293 PetscMPIInt commsize; 294 PetscBool localized,isascii; 295 PetscBool enable_mfem,enable_boundary,enable_ncmesh; 296 PetscBT pown,vown; 297 PetscErrorCode ierr; 298 PetscContainer glvis_container; 299 PetscBool periodic, enabled = PETSC_TRUE; 300 const char *fmt; 301 302 PetscFunctionBegin; 303 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 304 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 305 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 306 if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII"); 307 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr); 308 if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); 309 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 310 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 311 if (depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 312 313 /* get container: determines if a process visualizes is portion of the data or not */ 314 ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr); 315 if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); 316 { 317 PetscViewerGLVisInfo glvis_info; 318 ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr); 319 enabled = glvis_info->enabled; 320 fmt = glvis_info->fmt; 321 } 322 ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 323 ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 324 cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 325 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 326 ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr); 327 ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr); 328 if (periodic && !localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized"); 329 ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); 330 ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 331 ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr); 332 333 /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh 334 DMPlex does not currently support HO meshes, so there's no API for this */ 335 ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr); 336 337 /* 338 a couple of sections of the mesh specification are disabled 339 - boundary: unless we want to visualize boundary attributes or we have a periodic mesh 340 the boundary is not needed for proper mesh visualization 341 - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package 342 and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel) 343 */ 344 enable_boundary = periodic; 345 enable_ncmesh = PETSC_FALSE; 346 enable_mfem = PETSC_FALSE; 347 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr); 348 ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr); 349 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); 350 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); 351 ierr = PetscOptionsEnd();CHKERRQ(ierr); 352 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr); 353 if (enable_ncmesh && commsize > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel"); 354 355 /* Identify possible cells in the overlap */ 356 novl = 0; 357 pown = NULL; 358 if (commsize > 1) { 359 IS globalNum = NULL; 360 const PetscInt *gNum; 361 PetscBool ovl = PETSC_FALSE; 362 363 ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr); 364 ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 365 for (p=cStart; p<cEnd; p++) { 366 if (gNum[p-cStart] < 0) { 367 ovl = PETSC_TRUE; 368 novl++; 369 } 370 } 371 if (ovl) { 372 /* it may happen that pown get not destroyed, if the user closes the window while this function is running. 373 TODO: garbage collector? attach pown to dm? */ 374 ierr = PetscBTCreate(cEnd-cStart,&pown);CHKERRQ(ierr); 375 for (p=cStart; p<cEnd; p++) { 376 if (gNum[p-cStart] < 0) continue; 377 else { 378 ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr); 379 } 380 } 381 } 382 ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 383 } 384 385 /* return if this process is disabled */ 386 if (!enabled) { 387 ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 388 ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 389 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 390 ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 391 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 392 ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 393 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 394 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 395 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 396 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 397 ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 398 PetscFunctionReturn(0); 399 } 400 401 if (enable_mfem) { 402 if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */ 403 PetscInt vpc = 0; 404 char fec[64]; 405 int vids[8] = {0,1,2,3,4,5,6,7}; 406 int hexv[8] = {0,1,3,2,4,5,7,6},*dof; 407 PetscScalar *array,*ptr; 408 409 ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%dD_P1",dim);CHKERRQ(ierr); 410 if (cEnd-cStart) { 411 PetscInt fpc; 412 413 ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 414 switch(dim) { 415 case 1: 416 vpc = 2; 417 dof = hexv; 418 break; 419 case 2: 420 switch (fpc) { 421 case 3: 422 vpc = 3; 423 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 424 break; 425 case 4: 426 vpc = 4; 427 dof = hexv; 428 break; 429 default: 430 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 431 break; 432 } 433 break; 434 case 3: 435 switch (fpc) { 436 case 4: 437 vpc = 4; 438 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 439 break; 440 case 6: 441 vpc = 8; 442 dof = hexv; 443 break; 444 default: 445 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 446 break; 447 } 448 break; 449 default: 450 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 451 break; 452 } 453 ierr = DMPlexInvertCell(dim,vpc,vids);CHKERRQ(ierr); 454 } 455 ierr = VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);CHKERRQ(ierr); 456 ierr = PetscObjectCompose((PetscObject)dm,"_glvis_mesh_coords",(PetscObject)hovec);CHKERRQ(ierr); 457 ierr = PetscObjectDereference((PetscObject)hovec);CHKERRQ(ierr); 458 ierr = PetscObjectSetName((PetscObject)hovec,fec);CHKERRQ(ierr); 459 ierr = VecGetArray(hovec,&array);CHKERRQ(ierr); 460 ptr = array; 461 for (p=cStart;p<cEnd;p++) { 462 PetscInt csize,v,d; 463 PetscScalar *vals = NULL; 464 465 if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 466 ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 467 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); 468 for (v=0;v<vpc;v++) { 469 for (d=0;d<sdim;d++) { 470 ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d]; 471 } 472 } 473 ptr += vpc*sdim; 474 ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 475 } 476 ierr = VecRestoreArray(hovec,&array);CHKERRQ(ierr); 477 } 478 } 479 480 481 /* header */ 482 ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 483 484 /* topological dimension */ 485 ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 486 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 487 488 /* elements */ 489 /* TODO: is this the label we want to use for marking material IDs? 490 We should decide to have a single marker for these stuff 491 Proposal: DMSetMaterialLabel? 492 */ 493 ierr = DMGetLabel(dm,"Cell Sets",&label);CHKERRQ(ierr); 494 ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 495 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr); 496 for (p=cStart;p<cEnd;p++) { 497 int vids[8]; 498 PetscInt i,nv = 0,cid = -1,mid = 1; 499 500 if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 501 ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr); 502 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 503 ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 504 ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 505 for (i=0;i<nv;i++) { 506 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr); 507 } 508 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 509 } 510 511 /* boundary */ 512 ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 513 if (!enable_boundary) { 514 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 515 } else { 516 DMLabel perLabel; 517 PetscBT bfaces; 518 PetscInt fStart,fEnd,fEndInterior,*fcells; 519 PetscInt *faces = NULL,fpc = 0,vpf = 0, vpc = 0; 520 PetscInt fv1[] = {0,1}, 521 fv2tri[] = {0,1, 522 1,2, 523 2,0}, 524 fv2quad[] = {0,1, 525 1,2, 526 2,3, 527 3,0}, 528 fv3tet[] = {0,1,2, 529 0,3,1, 530 0,2,3, 531 2,1,3}, 532 fv3hex[] = {0,1,2,3, 533 4,5,6,7, 534 0,3,5,4, 535 2,1,7,6, 536 3,2,6,5, 537 0,4,7,1}; 538 539 /* determine orientation of boundary mesh */ 540 if (cEnd-cStart) { 541 ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 542 switch(dim) { 543 case 1: 544 if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc); 545 faces = fv1; 546 vpf = 1; 547 vpc = 2; 548 break; 549 case 2: 550 switch (fpc) { 551 case 3: 552 faces = fv2tri; 553 vpf = 2; 554 vpc = 3; 555 break; 556 case 4: 557 faces = fv2quad; 558 vpf = 2; 559 vpc = 4; 560 break; 561 default: 562 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 563 break; 564 } 565 break; 566 case 3: 567 switch (fpc) { 568 case 4: 569 faces = fv3tet; 570 vpf = 3; 571 vpc = 4; 572 break; 573 case 6: 574 faces = fv3hex; 575 vpf = 4; 576 vpc = 8; 577 break; 578 default: 579 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 580 break; 581 } 582 break; 583 default: 584 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 585 break; 586 } 587 } 588 ierr = DMPlexGetHybridBounds(dm,NULL,&fEndInterior,NULL,NULL);CHKERRQ(ierr); 589 ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr); 590 fEnd = fEndInterior < 0 ? fEnd : fEndInterior; 591 ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr); 592 ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr); 593 ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr); 594 ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 595 if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */ 596 ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr); 597 ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 598 ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr); 599 for (p=cStart;p<cEnd;p++) { 600 PetscInt dof; 601 ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr); 602 if (dof) { 603 PetscInt v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL; 604 PetscScalar *vals = NULL; 605 606 if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim); 607 if (dof/sdim != vpc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible number of cell dofs %D and vertices %D",dof/sdim,vpc); 608 ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 609 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 610 for (v=0;v<cellClosureSize;v++) 611 if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) { 612 vidxs = cellClosure + 2*v; 613 break; 614 } 615 if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices"); 616 for (v=0;v<vpc;v++) { 617 PetscInt s; 618 for (s=0;s<sdim;s++) { 619 if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+vpc*sdim])>PETSC_MACHINE_EPSILON) { 620 ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr); 621 } 622 } 623 } 624 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 625 ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 626 } 627 } 628 if (dim > 1) { 629 PetscInt eEnd,eStart,eEndInterior; 630 ierr = DMPlexGetHybridBounds(dm,NULL,NULL,&eEndInterior,NULL);CHKERRQ(ierr); 631 ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr); 632 eEnd = eEndInterior < 0 ? eEnd : eEndInterior; 633 for (p=eStart;p<eEnd;p++) { 634 const PetscInt *cone; 635 PetscInt coneSize,i; 636 PetscBool ispe = PETSC_TRUE; 637 638 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 639 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 640 for (i=0;i<coneSize;i++) { 641 PetscInt v; 642 643 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 644 ispe = (PetscBool)(ispe && (v==2)); 645 } 646 if (ispe && coneSize) { 647 ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 648 } 649 } 650 if (dim > 2) { 651 for (p=fStart;p<fEnd;p++) { 652 const PetscInt *cone; 653 PetscInt coneSize,i; 654 PetscBool ispe = PETSC_TRUE; 655 656 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 657 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 658 for (i=0;i<coneSize;i++) { 659 PetscInt v; 660 661 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 662 ispe = (PetscBool)(ispe && (v==2)); 663 } 664 if (ispe && coneSize) { 665 ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 666 } 667 } 668 } 669 } 670 } 671 for (p=fStart;p<fEnd;p++) { 672 const PetscInt *support; 673 PetscInt supportSize; 674 PetscBool isbf = PETSC_FALSE; 675 676 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 677 if (pown) { 678 PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE; 679 PetscInt i; 680 681 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 682 for (i=0;i<supportSize;i++) { 683 if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE; 684 else has_ghost = PETSC_TRUE; 685 } 686 isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)); 687 } else { 688 isbf = (PetscBool)(supportSize == 1); 689 } 690 if (!isbf && perLabel) { 691 const PetscInt *cone; 692 PetscInt coneSize,i; 693 694 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 695 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 696 isbf = PETSC_TRUE; 697 for (i=0;i<coneSize;i++) { 698 PetscInt v,d; 699 700 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 701 ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr); 702 isbf = (PetscBool)(isbf && v != d); 703 } 704 } 705 if (isbf) { 706 ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr); 707 } 708 } 709 /* count boundary faces */ 710 for (p=fStart,bf=0;p<fEnd;p++) { 711 if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 712 const PetscInt *support; 713 PetscInt supportSize,c; 714 715 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 716 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 717 if (pown) { 718 for (c=0;c<supportSize;c++) { 719 if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) { 720 bf++; 721 } 722 } 723 } else bf += supportSize; 724 } 725 } 726 /* TODO: is this the label we want to use for marking boundary facets? 727 We should decide to have a single marker for these stuff 728 Proposal: DMSetBoundaryLabel? 729 */ 730 ierr = DMGetLabel(dm,"marker",&label);CHKERRQ(ierr); 731 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr); 732 for (p=fStart;p<fEnd;p++) { 733 if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 734 const PetscInt *support; 735 PetscInt supportSize,c,nc = 0; 736 737 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 738 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 739 if (pown) { 740 for (c=0;c<supportSize;c++) { 741 if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) { 742 fcells[nc++] = support[c]; 743 } 744 } 745 } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c]; 746 for (c=0;c<nc;c++) { 747 const PetscInt *cone; 748 int vids[8]; 749 PetscInt i,cell,cl,nv,cid = -1,mid = -1; 750 751 cell = fcells[c]; 752 ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 753 for (cl=0;cl<fpc;cl++) 754 if (cone[cl] == p) 755 break; 756 757 /* face material id and type */ 758 ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr); 759 ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 760 /* vertex ids */ 761 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 762 for (i=0;i<vpf;i++) { 763 ierr = PetscViewerASCIIPrintf(viewer," %D",vids[faces[cl*vpf+i]]);CHKERRQ(ierr); 764 } 765 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 766 } 767 bf = bf-nc; 768 } 769 } 770 if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf); 771 ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr); 772 ierr = PetscFree(fcells);CHKERRQ(ierr); 773 } 774 775 /* mark owned vertices */ 776 vown = NULL; 777 if (pown) { 778 ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); 779 for (p=cStart;p<cEnd;p++) { 780 PetscInt i,closureSize,*closure = NULL; 781 782 if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue; 783 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 784 for (i=0;i<closureSize;i++) { 785 const PetscInt pp = closure[2*i]; 786 787 if (pp >= vStart && pp < vEnd) { 788 ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr); 789 } 790 } 791 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 792 } 793 } 794 795 /* vertex_parents (Non-conforming meshes) */ 796 parentSection = NULL; 797 if (enable_ncmesh) { 798 ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 799 } 800 if (parentSection) { 801 PetscInt vp,gvp; 802 803 for (vp=0,p=vStart;p<vEnd;p++) { 804 DMLabel dlabel; 805 PetscInt parent,depth; 806 807 if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 808 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 809 ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr); 810 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 811 if (parent != p) vp++; 812 } 813 ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 814 if (gvp) { 815 PetscInt maxsupp; 816 PetscBool *skip = NULL; 817 818 ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr); 819 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr); 820 ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr); 821 ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr); 822 for (p=vStart;p<vEnd;p++) { 823 DMLabel dlabel; 824 PetscInt parent; 825 826 if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 827 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 828 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 829 if (parent != p) { 830 int vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */ 831 PetscInt i,nv,size,n,numChildren,depth = -1; 832 const PetscInt *children; 833 ierr = DMPlexGetConeSize(dm,parent,&size);CHKERRQ(ierr); 834 switch (size) { 835 case 2: /* edge */ 836 nv = 0; 837 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 838 ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 839 ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr); 840 for (i=0;i<nv;i++) { 841 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr); 842 } 843 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 844 vp--; 845 break; 846 case 4: /* face */ 847 ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr); 848 for (n=0;n<numChildren;n++) { 849 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr); 850 if (!depth) { 851 const PetscInt *hvsupp,*hesupp,*cone; 852 PetscInt hvsuppSize,hesuppSize,coneSize; 853 PetscInt hv = children[n],he = -1,f; 854 855 ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr); 856 ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr); 857 ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr); 858 for (i=0;i<hvsuppSize;i++) { 859 PetscInt ep; 860 ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr); 861 if (ep != hvsupp[i]) { 862 he = hvsupp[i]; 863 } else { 864 skip[i] = PETSC_TRUE; 865 } 866 } 867 if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize); 868 ierr = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr); 869 vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 870 ierr = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr); 871 ierr = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr); 872 for (f=0;f<hesuppSize;f++) { 873 PetscInt j; 874 875 ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr); 876 ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr); 877 for (j=0;j<coneSize;j++) { 878 PetscInt k; 879 for (k=0;k<hvsuppSize;k++) { 880 if (hvsupp[k] == cone[j]) { 881 skip[k] = PETSC_TRUE; 882 break; 883 } 884 } 885 } 886 } 887 for (i=0;i<hvsuppSize;i++) { 888 if (!skip[i]) { 889 ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr); 890 vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 891 } 892 } 893 ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr); 894 for (i=0;i<2;i++) { 895 ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]-vStart);CHKERRQ(ierr); 896 } 897 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 898 vp--; 899 } 900 } 901 break; 902 default: 903 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %d",size); 904 } 905 } 906 } 907 ierr = PetscFree(skip);CHKERRQ(ierr); 908 } 909 if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp); 910 } 911 ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 912 ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); 913 914 /* vertices */ 915 if (hovec) { /* higher-order meshes */ 916 const char *fec; 917 PetscInt i,n; 918 919 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 920 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr); 921 ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr); 922 ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr); 923 ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr); 924 ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr); 925 ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr); 926 ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/ 927 ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr); 928 ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr); 929 if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim); 930 for (i=0;i<n/sdim;i++) { 931 PetscInt s; 932 for (s=0;s<sdim;s++) { 933 ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr); 934 } 935 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 936 } 937 ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr); 938 } else { 939 ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr); 940 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 941 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr); 942 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 943 ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr); 944 for (p=0;p<nvert/sdim;p++) { 945 PetscInt s; 946 for (s=0;s<sdim;s++) { 947 ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[p*sdim+s]));CHKERRQ(ierr); 948 } 949 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 950 } 951 ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr); 952 } 953 PetscFunctionReturn(0); 954 } 955 956 /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */ 957 PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer) 958 { 959 PetscErrorCode ierr; 960 PetscBool isglvis,isascii; 961 962 PetscFunctionBegin; 963 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 964 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 965 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr); 966 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 967 if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII"); 968 if (isglvis) { 969 PetscViewer view; 970 PetscViewerGLVisType type; 971 972 ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr); 973 ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr); 974 if (view) { /* in the socket case, it may happen that the connection failed */ 975 if (type == PETSC_VIEWER_GLVIS_SOCKET) { 976 PetscMPIInt size,rank; 977 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 978 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 979 ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr); 980 } 981 ierr = DMPlexView_GLVis_ASCII(dm,view);CHKERRQ(ierr); 982 ierr = PetscViewerFlush(view);CHKERRQ(ierr); 983 if (type == PETSC_VIEWER_GLVIS_SOCKET) { 984 PetscInt dim; 985 const char* name; 986 987 ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr); 988 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 989 ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr); 990 ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr); 991 } 992 } 993 } else { 994 ierr = DMPlexView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr); 995 } 996 PetscFunctionReturn(0); 997 } 998