18135c375SStefano Zampini #include <petsc/private/glvisviewerimpl.h> 28135c375SStefano Zampini #include <petsc/private/petscimpl.h> 38135c375SStefano Zampini #include <petsc/private/dmpleximpl.h> 48135c375SStefano Zampini #include <petscbt.h> 58135c375SStefano Zampini #include <petscdmplex.h> 68135c375SStefano Zampini #include <petscsf.h> 78135c375SStefano Zampini #include <petscds.h> 88135c375SStefano Zampini 98135c375SStefano Zampini typedef struct { 108135c375SStefano Zampini PetscInt nf; 118135c375SStefano Zampini VecScatter *scctx; 128135c375SStefano Zampini } GLVisViewerCtx; 138135c375SStefano Zampini 148135c375SStefano Zampini static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) 158135c375SStefano Zampini { 168135c375SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 178135c375SStefano Zampini PetscInt i; 188135c375SStefano Zampini PetscErrorCode ierr; 198135c375SStefano Zampini 208135c375SStefano Zampini PetscFunctionBegin; 218135c375SStefano Zampini for (i=0;i<ctx->nf;i++) { 228135c375SStefano Zampini ierr = VecScatterDestroy(&ctx->scctx[i]);CHKERRQ(ierr); 238135c375SStefano Zampini } 248135c375SStefano Zampini ierr = PetscFree(ctx->scctx);CHKERRQ(ierr); 2551f51421SSatish Balay ierr = PetscFree(vctx);CHKERRQ(ierr); 268135c375SStefano Zampini PetscFunctionReturn(0); 278135c375SStefano Zampini } 288135c375SStefano Zampini 298135c375SStefano Zampini static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) 308135c375SStefano Zampini { 318135c375SStefano Zampini GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; 328135c375SStefano Zampini PetscInt f; 338135c375SStefano Zampini PetscErrorCode ierr; 348135c375SStefano Zampini 358135c375SStefano Zampini PetscFunctionBegin; 368135c375SStefano Zampini for (f=0;f<nf;f++) { 378135c375SStefano Zampini ierr = VecScatterBegin(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 388135c375SStefano Zampini ierr = VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 398135c375SStefano Zampini } 408135c375SStefano Zampini PetscFunctionReturn(0); 418135c375SStefano Zampini } 428135c375SStefano Zampini 438135c375SStefano Zampini /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */ 448135c375SStefano Zampini PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer) 458135c375SStefano Zampini { 468135c375SStefano Zampini DM dm = (DM)odm; 478135c375SStefano Zampini Vec xlocal,xfield; 488135c375SStefano Zampini PetscDS ds; 498135c375SStefano Zampini IS globalNum,isfield; 508135c375SStefano Zampini PetscBT vown; 518135c375SStefano Zampini char **fieldname = NULL,**fec_type = NULL; 528135c375SStefano Zampini const PetscInt *gNum; 53bb77a09fSStefano Zampini PetscInt *nlocal,*bs,*idxs,*dims; 548135c375SStefano Zampini PetscInt f,maxfields,nfields,c,totc,totdofs,Nv,cum,i; 5511a4995dSStefano Zampini PetscInt dim,cStart,cEnd,cEndInterior,vStart,vEnd; 568135c375SStefano Zampini GLVisViewerCtx *ctx; 578135c375SStefano Zampini PetscSection s; 588135c375SStefano Zampini PetscErrorCode ierr; 598135c375SStefano Zampini 608135c375SStefano Zampini PetscFunctionBegin; 618135c375SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 628135c375SStefano Zampini ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 638135c375SStefano Zampini ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 648135c375SStefano Zampini ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 658135c375SStefano Zampini cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 668135c375SStefano Zampini ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr); 678135c375SStefano Zampini ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 688135c375SStefano Zampini ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); 698135c375SStefano Zampini for (c = cStart, totc = 0; c < cEnd; c++) { 708135c375SStefano Zampini if (gNum[c-cStart] >= 0) { 718135c375SStefano Zampini PetscInt i,numPoints,*points = NULL; 728135c375SStefano Zampini 738135c375SStefano Zampini totc++; 748135c375SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 758135c375SStefano Zampini for (i=0;i<numPoints*2;i+= 2) { 768135c375SStefano Zampini if ((points[i] >= vStart) && (points[i] < vEnd)) { 778135c375SStefano Zampini ierr = PetscBTSet(vown,points[i]-vStart);CHKERRQ(ierr); 788135c375SStefano Zampini } 798135c375SStefano Zampini } 808135c375SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 818135c375SStefano Zampini } 828135c375SStefano Zampini } 8377eacf09SStefano Zampini for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscLikely(PetscBTLookup(vown,f))) Nv++; 848135c375SStefano Zampini 858135c375SStefano Zampini ierr = DMCreateLocalVector(dm,&xlocal);CHKERRQ(ierr); 868135c375SStefano Zampini ierr = VecGetLocalSize(xlocal,&totdofs);CHKERRQ(ierr); 878135c375SStefano Zampini ierr = DMGetDefaultSection(dm,&s);CHKERRQ(ierr); 888135c375SStefano Zampini ierr = PetscSectionGetNumFields(s,&nfields);CHKERRQ(ierr); 898135c375SStefano Zampini for (f=0,maxfields=0;f<nfields;f++) { 908135c375SStefano Zampini PetscInt bs; 918135c375SStefano Zampini 928135c375SStefano Zampini ierr = PetscSectionGetFieldComponents(s,f,&bs);CHKERRQ(ierr); 938135c375SStefano Zampini maxfields += bs; 948135c375SStefano Zampini } 95bb77a09fSStefano Zampini ierr = PetscCalloc6(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs);CHKERRQ(ierr); 968135c375SStefano Zampini ierr = PetscNew(&ctx);CHKERRQ(ierr); 978135c375SStefano Zampini ierr = PetscCalloc1(maxfields,&ctx->scctx);CHKERRQ(ierr); 988135c375SStefano Zampini ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 998135c375SStefano Zampini if (ds) { 1008135c375SStefano Zampini for (f=0;f<nfields;f++) { 1018135c375SStefano Zampini const char* fname; 1028135c375SStefano Zampini char name[256]; 1038135c375SStefano Zampini PetscObject disc; 1048135c375SStefano Zampini size_t len; 1058135c375SStefano Zampini 1068135c375SStefano Zampini ierr = PetscSectionGetFieldName(s,f,&fname);CHKERRQ(ierr); 1078135c375SStefano Zampini ierr = PetscStrlen(fname,&len);CHKERRQ(ierr); 1088135c375SStefano Zampini if (len) { 1098135c375SStefano Zampini ierr = PetscStrcpy(name,fname);CHKERRQ(ierr); 1108135c375SStefano Zampini } else { 1118135c375SStefano Zampini ierr = PetscSNPrintf(name,256,"Field%d",f);CHKERRQ(ierr); 1128135c375SStefano Zampini } 1138135c375SStefano Zampini ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr); 1148135c375SStefano Zampini if (disc) { 1158135c375SStefano Zampini PetscClassId id; 1168135c375SStefano Zampini PetscInt Nc; 1178135c375SStefano Zampini char fec[64]; 1188135c375SStefano Zampini 1198135c375SStefano Zampini ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr); 1208135c375SStefano Zampini if (id == PETSCFE_CLASSID) { 1218135c375SStefano Zampini PetscFE fem = (PetscFE)disc; 1228135c375SStefano Zampini PetscDualSpace sp; 1238135c375SStefano Zampini PetscDualSpaceType spname; 1248135c375SStefano Zampini PetscInt order; 1258135c375SStefano Zampini PetscBool islag,continuous,H1 = PETSC_TRUE; 1268135c375SStefano Zampini 1278135c375SStefano Zampini ierr = PetscFEGetNumComponents(fem,&Nc);CHKERRQ(ierr); 1288135c375SStefano Zampini ierr = PetscFEGetDualSpace(fem,&sp);CHKERRQ(ierr); 1298135c375SStefano Zampini ierr = PetscDualSpaceGetType(sp,&spname);CHKERRQ(ierr); 1308135c375SStefano Zampini ierr = PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);CHKERRQ(ierr); 1318135c375SStefano Zampini if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space"); 1328135c375SStefano Zampini ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr); 1338135c375SStefano Zampini if (!continuous) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported"); 1348135c375SStefano Zampini ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 1358135c375SStefano Zampini if (continuous && order > 0) { 1368135c375SStefano Zampini ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%dD_P1",dim);CHKERRQ(ierr); 1378135c375SStefano Zampini } else { 1388135c375SStefano Zampini H1 = PETSC_FALSE; 1398135c375SStefano Zampini ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P%d",dim,order);CHKERRQ(ierr); 1408135c375SStefano Zampini } 1418135c375SStefano Zampini ierr = PetscStrallocpy(name,&fieldname[ctx->nf]);CHKERRQ(ierr); 1428135c375SStefano Zampini bs[ctx->nf] = Nc; 143bb77a09fSStefano Zampini dims[ctx->nf] = dim; 1448135c375SStefano Zampini if (H1) { 1458135c375SStefano Zampini nlocal[ctx->nf] = Nc * Nv; 1468135c375SStefano Zampini ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 1478135c375SStefano Zampini ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr); 1488135c375SStefano Zampini for (i=0,cum=0;i<vEnd-vStart;i++) { 1498135c375SStefano Zampini PetscInt j,off; 1508135c375SStefano Zampini 1518135c375SStefano Zampini if (PetscUnlikely(!PetscBTLookup(vown,i))) continue; 1528135c375SStefano Zampini ierr = PetscSectionGetFieldOffset(s,i+vStart,f,&off);CHKERRQ(ierr); 1538135c375SStefano Zampini for (j=0;j<Nc;j++) idxs[cum++] = off + j; 1548135c375SStefano Zampini } 1558135c375SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 1568135c375SStefano Zampini } else { 1578135c375SStefano Zampini nlocal[ctx->nf] = Nc * totc; 1588135c375SStefano Zampini ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 1598135c375SStefano Zampini ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr); 1608135c375SStefano Zampini for (i=0,cum=0;i<cEnd-cStart;i++) { 1618135c375SStefano Zampini PetscInt j,off; 1628135c375SStefano Zampini 1638135c375SStefano Zampini if (PetscUnlikely(gNum[i] < 0)) continue; 1648135c375SStefano Zampini ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr); 1658135c375SStefano Zampini for (j=0;j<Nc;j++) idxs[cum++] = off + j; 1668135c375SStefano Zampini } 1678135c375SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 1688135c375SStefano Zampini } 1698135c375SStefano Zampini ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr); 1708135c375SStefano Zampini ierr = VecDestroy(&xfield);CHKERRQ(ierr); 1718135c375SStefano Zampini ierr = ISDestroy(&isfield);CHKERRQ(ierr); 1728135c375SStefano Zampini ctx->nf++; 1738135c375SStefano Zampini } else if (id == PETSCFV_CLASSID) { 1748135c375SStefano Zampini PetscInt c; 1758135c375SStefano Zampini 1768135c375SStefano Zampini ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr); 1778135c375SStefano Zampini ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%dD_P0",dim);CHKERRQ(ierr); 1788135c375SStefano Zampini for (c = 0; c < Nc; c++) { 1798135c375SStefano Zampini char comp[256]; 1808135c375SStefano Zampini ierr = PetscSNPrintf(comp,256,"%s-Comp%d",name,c);CHKERRQ(ierr); 1818135c375SStefano Zampini ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr); 1828135c375SStefano Zampini bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */ 1838135c375SStefano Zampini nlocal[ctx->nf] = totc; 184bb77a09fSStefano Zampini dims[ctx->nf] = dim; 1858135c375SStefano Zampini ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 1868135c375SStefano Zampini ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr); 1878135c375SStefano Zampini for (i=0,cum=0;i<cEnd-cStart;i++) { 1888135c375SStefano Zampini PetscInt off; 1898135c375SStefano Zampini 1908135c375SStefano Zampini if (PetscUnlikely(gNum[i])<0) continue; 1918135c375SStefano Zampini ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr); 1928135c375SStefano Zampini idxs[cum++] = off + c; 1938135c375SStefano Zampini } 1948135c375SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 1958135c375SStefano Zampini ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr); 1968135c375SStefano Zampini ierr = VecDestroy(&xfield);CHKERRQ(ierr); 1978135c375SStefano Zampini ierr = ISDestroy(&isfield);CHKERRQ(ierr); 1988135c375SStefano Zampini ctx->nf++; 1998135c375SStefano Zampini } 2008135c375SStefano Zampini } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f); 2018135c375SStefano Zampini } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f); 2028135c375SStefano Zampini } 2038135c375SStefano Zampini } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM"); 2048135c375SStefano Zampini ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); 2058135c375SStefano Zampini ierr = VecDestroy(&xlocal);CHKERRQ(ierr); 2068135c375SStefano Zampini ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 2078135c375SStefano Zampini 2088135c375SStefano Zampini /* customize the viewer */ 209bb77a09fSStefano Zampini ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fieldname,(const char**)fec_type,nlocal,bs,dims,DMPlexSampleGLVisFields_Private,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr); 2108135c375SStefano Zampini for (f=0;f<ctx->nf;f++) { 2118135c375SStefano Zampini ierr = PetscFree(fieldname[f]);CHKERRQ(ierr); 2128135c375SStefano Zampini ierr = PetscFree(fec_type[f]);CHKERRQ(ierr); 2138135c375SStefano Zampini } 214bb77a09fSStefano Zampini ierr = PetscFree6(fieldname,nlocal,bs,dims,fec_type,idxs);CHKERRQ(ierr); 2158135c375SStefano Zampini PetscFunctionReturn(0); 2168135c375SStefano Zampini } 2178135c375SStefano Zampini 2188135c375SStefano Zampini typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_UNDEF} MFEM_cid; 2198135c375SStefano Zampini 2208135c375SStefano Zampini MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF, MFEM_UNDEF}, 2218135c375SStefano Zampini {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF, MFEM_UNDEF}, 2228135c375SStefano Zampini {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF, MFEM_UNDEF}, 2238135c375SStefano Zampini {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_UNDEF, MFEM_CUBE } }; 2248135c375SStefano Zampini 225cc0d3ed7SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt p, PetscInt *mid, PetscInt *cid) 2268135c375SStefano Zampini { 2278135c375SStefano Zampini DMLabel dlabel; 2288135c375SStefano Zampini PetscInt depth,csize; 2298135c375SStefano Zampini PetscErrorCode ierr; 2308135c375SStefano Zampini 2318135c375SStefano Zampini PetscFunctionBegin; 2328135c375SStefano Zampini ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 2338135c375SStefano Zampini ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr); 2348135c375SStefano Zampini ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr); 2358135c375SStefano Zampini if (label) { 2368135c375SStefano Zampini ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr); 23777eacf09SStefano Zampini } else *mid = 1; 2388135c375SStefano Zampini *cid = mfem_table_cid[depth][csize]; 2398135c375SStefano Zampini PetscFunctionReturn(0); 2408135c375SStefano Zampini } 2418135c375SStefano Zampini 242cc0d3ed7SStefano Zampini static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, int vids[]) 2438135c375SStefano Zampini { 244cc0d3ed7SStefano Zampini PetscInt dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL; 2458135c375SStefano Zampini PetscErrorCode ierr; 2468135c375SStefano Zampini 2478135c375SStefano Zampini PetscFunctionBegin; 2488135c375SStefano Zampini ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 249cc0d3ed7SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 250cc0d3ed7SStefano Zampini sdim = dim; 251cc0d3ed7SStefano Zampini if (csec) { 252cc0d3ed7SStefano Zampini ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 253cc0d3ed7SStefano Zampini ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr); 254cc0d3ed7SStefano Zampini off = off/sdim; 255cc0d3ed7SStefano Zampini ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 256cc0d3ed7SStefano Zampini } 257cc0d3ed7SStefano Zampini if (!dof) { 2588135c375SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 2598135c375SStefano Zampini for (i=0,q=0;i<numPoints*2;i+= 2) 2608135c375SStefano Zampini if ((points[i] >= vStart) && (points[i] < vEnd)) 261cc0d3ed7SStefano Zampini vids[q++] = points[i]-vStart+off; 2628135c375SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 263cc0d3ed7SStefano Zampini } else { 264cc0d3ed7SStefano Zampini ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr); 265cc0d3ed7SStefano Zampini ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 266cc0d3ed7SStefano Zampini for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q; 267cc0d3ed7SStefano Zampini } 2688135c375SStefano Zampini *nv = q; 2698135c375SStefano Zampini PetscFunctionReturn(0); 2708135c375SStefano Zampini } 2718135c375SStefano Zampini 27277eacf09SStefano Zampini /* 27377eacf09SStefano Zampini ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR 27477eacf09SStefano Zampini Higher order meshes are also supported 27577eacf09SStefano Zampini */ 2768135c375SStefano Zampini static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer) 2778135c375SStefano Zampini { 2788135c375SStefano Zampini DMLabel label; 2798135c375SStefano Zampini PetscSection coordSection,parentSection; 28077eacf09SStefano Zampini Vec coordinates,hovec; 2818135c375SStefano Zampini IS globalNum = NULL; 2828135c375SStefano Zampini const PetscScalar *array; 2838135c375SStefano Zampini const PetscInt *gNum; 2848135c375SStefano Zampini PetscInt bf,p,sdim,dim,depth,novl; 285cc0d3ed7SStefano Zampini PetscInt cStart,cEnd,cEndInterior,vStart,vEnd,nvert; 2868135c375SStefano Zampini PetscMPIInt commsize; 287cc0d3ed7SStefano Zampini PetscBool localized,ovl,isascii,enable_boundary,enable_ncmesh; 2888135c375SStefano Zampini PetscBT pown; 2898135c375SStefano Zampini PetscErrorCode ierr; 2908135c375SStefano Zampini PetscContainer glvis_container; 29177eacf09SStefano Zampini PetscBool periodic, enabled = PETSC_TRUE; 29277eacf09SStefano Zampini const char *fmt; 2938135c375SStefano Zampini 2948135c375SStefano Zampini PetscFunctionBegin; 2958135c375SStefano Zampini PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2968135c375SStefano Zampini PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2978135c375SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 2988135c375SStefano Zampini if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII"); 2998135c375SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr); 3008135c375SStefano Zampini if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); 3018135c375SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 3028135c375SStefano Zampini ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 3038135c375SStefano Zampini if (depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated"); 3048135c375SStefano Zampini 3058135c375SStefano Zampini /* get container: determines if a process visualizes is portion of the data or not */ 3068135c375SStefano Zampini ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr); 3078135c375SStefano Zampini if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); 3088135c375SStefano Zampini { 3098135c375SStefano Zampini PetscViewerGLVisInfo glvis_info; 3108135c375SStefano Zampini ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr); 3118135c375SStefano Zampini enabled = glvis_info->enabled; 31277eacf09SStefano Zampini fmt = glvis_info->fmt; 3138135c375SStefano Zampini } 314cc0d3ed7SStefano Zampini ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 3158135c375SStefano Zampini ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 3168135c375SStefano Zampini cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 3178135c375SStefano Zampini ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 31877eacf09SStefano Zampini ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr); 319cc0d3ed7SStefano Zampini ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr); 32077eacf09SStefano Zampini if (periodic && !localized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized"); 321cc0d3ed7SStefano Zampini ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); 32277eacf09SStefano Zampini ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 32377eacf09SStefano Zampini ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr); 32477eacf09SStefano Zampini 32577eacf09SStefano Zampini /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh 32677eacf09SStefano Zampini DMPlex does not currently support HO meshes, so there's no API for this */ 32777eacf09SStefano Zampini ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr); 3288135c375SStefano Zampini 3298135c375SStefano Zampini /* 3308135c375SStefano Zampini a couple of sections of the mesh specification are disabled 33177eacf09SStefano Zampini - boundary: unless we want to visualize boundary attributes or we have a periodic mesh 33277eacf09SStefano Zampini the boundary is not needed for proper mesh visualization 33377eacf09SStefano Zampini - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package 33477eacf09SStefano Zampini and be able to derefine the mesh 3358135c375SStefano Zampini */ 336*e74666bcSStefano Zampini enable_boundary = periodic; 3378135c375SStefano Zampini enable_ncmesh = PETSC_FALSE; 3388135c375SStefano Zampini ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr); 3398135c375SStefano Zampini ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation; useful for debugging purposes",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr); 3408135c375SStefano Zampini ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation; useful for debugging purposes",NULL,enable_ncmesh,&enable_ncmesh,NULL);CHKERRQ(ierr); 3418135c375SStefano Zampini ierr = PetscOptionsEnd();CHKERRQ(ierr); 3428135c375SStefano Zampini 3438135c375SStefano Zampini /* Identify possible cells in the overlap */ 3448135c375SStefano Zampini gNum = NULL; 3458135c375SStefano Zampini novl = 0; 3468135c375SStefano Zampini ovl = PETSC_FALSE; 3478135c375SStefano Zampini pown = NULL; 3488135c375SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr); 3498135c375SStefano Zampini if (commsize > 1) { 3508135c375SStefano Zampini ierr = DMPlexGetCellNumbering(dm,&globalNum);CHKERRQ(ierr); 3518135c375SStefano Zampini ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 3528135c375SStefano Zampini for (p=cStart; p<cEnd; p++) { 3538135c375SStefano Zampini if (gNum[p-cStart] < 0) { 3548135c375SStefano Zampini ovl = PETSC_TRUE; 3558135c375SStefano Zampini novl++; 3568135c375SStefano Zampini } 3578135c375SStefano Zampini } 3588135c375SStefano Zampini if (ovl) { 3598135c375SStefano Zampini /* it may happen that pown get not destroyed, if the user closes the window while this function is running. 3608135c375SStefano Zampini TODO: garbage collector? attach pown to dm? */ 3618135c375SStefano Zampini ierr = PetscBTCreate(PetscMax(cEnd-cStart,vEnd-vStart),&pown);CHKERRQ(ierr); 3628135c375SStefano Zampini } 3638135c375SStefano Zampini } 3648135c375SStefano Zampini 3658135c375SStefano Zampini if (!enabled) { 3668135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 3678135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 3688135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 3698135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 3708135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 3718135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 3728135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 3738135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 3748135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 3758135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 3768135c375SStefano Zampini ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 3778135c375SStefano Zampini if (globalNum) { 3788135c375SStefano Zampini ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 3798135c375SStefano Zampini } 3808135c375SStefano Zampini PetscFunctionReturn(0); 3818135c375SStefano Zampini } 3828135c375SStefano Zampini 3838135c375SStefano Zampini /* header */ 3848135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 3858135c375SStefano Zampini 3868135c375SStefano Zampini /* topological dimension */ 3878135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 3888135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 3898135c375SStefano Zampini 3908135c375SStefano Zampini /* elements */ 3918135c375SStefano Zampini /* TODO: is this the label we want to use for marking material IDs? 3928135c375SStefano Zampini We should decide to have a single marker for these stuff 3938135c375SStefano Zampini Proposal: DMSetMaterialLabel? 3948135c375SStefano Zampini */ 395d3e05d38SLisandro Dalcin ierr = DMGetLabel(dm,"Cell Sets",&label);CHKERRQ(ierr); 3968135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 3978135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr); 3988135c375SStefano Zampini for (p=cStart;p<cEnd;p++) { 39911a4995dSStefano Zampini int vids[8]; 40011a4995dSStefano Zampini PetscInt i,nv = 0,cid = -1,mid = 1; 4018135c375SStefano Zampini 4028135c375SStefano Zampini if (ovl) { 4038135c375SStefano Zampini if (gNum[p-cStart] < 0) continue; 4048135c375SStefano Zampini else { 4058135c375SStefano Zampini ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr); 4068135c375SStefano Zampini } 4078135c375SStefano Zampini } 4088135c375SStefano Zampini ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr); 40977eacf09SStefano Zampini ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 41077eacf09SStefano Zampini ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 4118135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 4128135c375SStefano Zampini for (i=0;i<nv;i++) { 41311a4995dSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr); 4148135c375SStefano Zampini } 4158135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 4168135c375SStefano Zampini } 4178135c375SStefano Zampini 418cc0d3ed7SStefano Zampini /* boundary */ 419cc0d3ed7SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 420cc0d3ed7SStefano Zampini if (!enable_boundary) { 421cc0d3ed7SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 422cc0d3ed7SStefano Zampini } else { 42377eacf09SStefano Zampini DMLabel perLabel; 42477eacf09SStefano Zampini PetscBT bfaces; 42577eacf09SStefano Zampini PetscInt fStart,fEnd,fEndInterior,*fcells; 42677eacf09SStefano Zampini PetscInt *faces = NULL,fpc = 0,vpf = 0, vpc = 0; 427cc0d3ed7SStefano Zampini PetscInt fv1[] = {0,1}, 428cc0d3ed7SStefano Zampini fv2tri[] = {0,1, 429cc0d3ed7SStefano Zampini 1,2, 430cc0d3ed7SStefano Zampini 2,0}, 431cc0d3ed7SStefano Zampini fv2quad[] = {0,1, 432cc0d3ed7SStefano Zampini 1,2, 433cc0d3ed7SStefano Zampini 2,3, 434cc0d3ed7SStefano Zampini 3,0}, 435cc0d3ed7SStefano Zampini fv3tet[] = {0,1,2, 436cc0d3ed7SStefano Zampini 0,3,1, 437cc0d3ed7SStefano Zampini 0,2,3, 438cc0d3ed7SStefano Zampini 2,1,3}, 439cc0d3ed7SStefano Zampini fv3hex[] = {0,1,2,3, 440cc0d3ed7SStefano Zampini 4,5,6,7, 441cc0d3ed7SStefano Zampini 0,3,5,4, 442cc0d3ed7SStefano Zampini 2,1,7,6, 443cc0d3ed7SStefano Zampini 3,2,6,5, 444cc0d3ed7SStefano Zampini 0,4,7,1}; 445cc0d3ed7SStefano Zampini 4468135c375SStefano Zampini /* determine orientation of boundary mesh */ 4478135c375SStefano Zampini if (cEnd-cStart) { 4488135c375SStefano Zampini ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 4498135c375SStefano Zampini switch(dim) { 4508135c375SStefano Zampini case 1: 4518135c375SStefano Zampini if (fpc != 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case faces per cell %D",fpc); 4528135c375SStefano Zampini faces = fv1; 4538135c375SStefano Zampini vpf = 1; 45477eacf09SStefano Zampini vpc = 2; 4558135c375SStefano Zampini break; 4568135c375SStefano Zampini case 2: 4578135c375SStefano Zampini switch (fpc) { 4588135c375SStefano Zampini case 3: 4598135c375SStefano Zampini faces = fv2tri; 4608135c375SStefano Zampini vpf = 2; 46177eacf09SStefano Zampini vpc = 3; 4628135c375SStefano Zampini break; 4638135c375SStefano Zampini case 4: 4648135c375SStefano Zampini faces = fv2quad; 4658135c375SStefano Zampini vpf = 2; 46677eacf09SStefano Zampini vpc = 4; 4678135c375SStefano Zampini break; 4688135c375SStefano Zampini default: 4698135c375SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 4708135c375SStefano Zampini break; 4718135c375SStefano Zampini } 4728135c375SStefano Zampini break; 4738135c375SStefano Zampini case 3: 4748135c375SStefano Zampini switch (fpc) { 4758135c375SStefano Zampini case 4: 4768135c375SStefano Zampini faces = fv3tet; 4778135c375SStefano Zampini vpf = 3; 47877eacf09SStefano Zampini vpc = 4; 4798135c375SStefano Zampini break; 4808135c375SStefano Zampini case 6: 4818135c375SStefano Zampini faces = fv3hex; 4828135c375SStefano Zampini vpf = 4; 48377eacf09SStefano Zampini vpc = 8; 4848135c375SStefano Zampini break; 4858135c375SStefano Zampini default: 4868135c375SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 4878135c375SStefano Zampini break; 4888135c375SStefano Zampini } 4898135c375SStefano Zampini break; 4908135c375SStefano Zampini default: 4918135c375SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 4928135c375SStefano Zampini break; 4938135c375SStefano Zampini } 4948135c375SStefano Zampini } 495cc0d3ed7SStefano Zampini ierr = DMPlexGetHybridBounds(dm,NULL,&fEndInterior,NULL,NULL);CHKERRQ(ierr); 4968135c375SStefano Zampini ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr); 4978135c375SStefano Zampini fEnd = fEndInterior < 0 ? fEnd : fEndInterior; 49877eacf09SStefano Zampini ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr); 49977eacf09SStefano Zampini ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr); 50077eacf09SStefano Zampini ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr); 50177eacf09SStefano Zampini ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 50277eacf09SStefano Zampini if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */ 50377eacf09SStefano Zampini ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr); 50477eacf09SStefano Zampini ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 50577eacf09SStefano Zampini ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr); 50677eacf09SStefano Zampini for (p=cStart;p<cEnd;p++) { 50777eacf09SStefano Zampini PetscInt dof; 50877eacf09SStefano Zampini ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr); 50977eacf09SStefano Zampini if (dof) { 51077eacf09SStefano Zampini PetscInt v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL; 51177eacf09SStefano Zampini PetscScalar *vals = NULL; 51277eacf09SStefano Zampini 51377eacf09SStefano Zampini if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim); 51477eacf09SStefano Zampini if (dof/sdim != vpc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible number of cell dofs %D and vertices %D",dof/sdim,vpc); 51577eacf09SStefano Zampini ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 51677eacf09SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 51777eacf09SStefano Zampini for (v=0;v<cellClosureSize;v++) 51877eacf09SStefano Zampini if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) { 51977eacf09SStefano Zampini vidxs = cellClosure + 2*v; 52077eacf09SStefano Zampini break; 52177eacf09SStefano Zampini } 52277eacf09SStefano Zampini if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices"); 52377eacf09SStefano Zampini for (v=0;v<vpc;v++) { 52477eacf09SStefano Zampini PetscInt s; 52577eacf09SStefano Zampini for (s=0;s<sdim;s++) { 52677eacf09SStefano Zampini if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+vpc*sdim])>PETSC_MACHINE_EPSILON) { 52777eacf09SStefano Zampini ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr); 52877eacf09SStefano Zampini } 52977eacf09SStefano Zampini } 53077eacf09SStefano Zampini } 53177eacf09SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 53277eacf09SStefano Zampini ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 53377eacf09SStefano Zampini } 53477eacf09SStefano Zampini } 53577eacf09SStefano Zampini if (dim > 1) { 53677eacf09SStefano Zampini PetscInt eEnd,eStart,eEndInterior; 53777eacf09SStefano Zampini ierr = DMPlexGetHybridBounds(dm,NULL,NULL,&eEndInterior,NULL);CHKERRQ(ierr); 53877eacf09SStefano Zampini ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr); 53977eacf09SStefano Zampini eEnd = eEndInterior < 0 ? eEnd : eEndInterior; 54077eacf09SStefano Zampini for (p=eStart;p<eEnd;p++) { 54177eacf09SStefano Zampini const PetscInt *cone; 54277eacf09SStefano Zampini PetscInt coneSize,i; 54377eacf09SStefano Zampini PetscBool ispe = PETSC_TRUE; 54477eacf09SStefano Zampini 54577eacf09SStefano Zampini ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 54677eacf09SStefano Zampini ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 54777eacf09SStefano Zampini for (i=0;i<coneSize;i++) { 54877eacf09SStefano Zampini PetscInt v; 54977eacf09SStefano Zampini 55077eacf09SStefano Zampini ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 55177eacf09SStefano Zampini ispe = (PetscBool)(ispe && (v==2)); 55277eacf09SStefano Zampini } 55377eacf09SStefano Zampini if (ispe && coneSize) { 55477eacf09SStefano Zampini ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 55577eacf09SStefano Zampini } 55677eacf09SStefano Zampini } 55777eacf09SStefano Zampini if (dim > 2) { 55877eacf09SStefano Zampini for (p=fStart;p<fEnd;p++) { 55977eacf09SStefano Zampini const PetscInt *cone; 56077eacf09SStefano Zampini PetscInt coneSize,i; 56177eacf09SStefano Zampini PetscBool ispe = PETSC_TRUE; 56277eacf09SStefano Zampini 56377eacf09SStefano Zampini ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 56477eacf09SStefano Zampini ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 56577eacf09SStefano Zampini for (i=0;i<coneSize;i++) { 56677eacf09SStefano Zampini PetscInt v; 56777eacf09SStefano Zampini 56877eacf09SStefano Zampini ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 56977eacf09SStefano Zampini ispe = (PetscBool)(ispe && (v==2)); 57077eacf09SStefano Zampini } 57177eacf09SStefano Zampini if (ispe && coneSize) { 57277eacf09SStefano Zampini ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 57377eacf09SStefano Zampini } 57477eacf09SStefano Zampini } 57577eacf09SStefano Zampini } 57677eacf09SStefano Zampini } 57777eacf09SStefano Zampini } 57877eacf09SStefano Zampini for (p=fStart;p<fEnd;p++) { 57977eacf09SStefano Zampini const PetscInt *support; 5808135c375SStefano Zampini PetscInt supportSize; 58177eacf09SStefano Zampini PetscBool isbf = PETSC_FALSE; 5828135c375SStefano Zampini 5838135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 58477eacf09SStefano Zampini if (ovl) { 5858135c375SStefano Zampini PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE; 58677eacf09SStefano Zampini PetscInt i; 58777eacf09SStefano Zampini 58877eacf09SStefano Zampini ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 58977eacf09SStefano Zampini for (i=0;i<supportSize;i++) { 59077eacf09SStefano Zampini if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE; 59177eacf09SStefano Zampini else has_ghost = PETSC_TRUE; 59277eacf09SStefano Zampini } 59377eacf09SStefano Zampini isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)); 59477eacf09SStefano Zampini } else { 59577eacf09SStefano Zampini isbf = (PetscBool)(supportSize == 1); 59677eacf09SStefano Zampini } 59777eacf09SStefano Zampini if (!isbf && perLabel) { 59877eacf09SStefano Zampini const PetscInt *cone; 59977eacf09SStefano Zampini PetscInt coneSize,i; 60077eacf09SStefano Zampini 60177eacf09SStefano Zampini ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 60277eacf09SStefano Zampini ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 60377eacf09SStefano Zampini isbf = PETSC_TRUE; 60477eacf09SStefano Zampini for (i=0;i<coneSize;i++) { 60577eacf09SStefano Zampini PetscInt v,d; 60677eacf09SStefano Zampini 60777eacf09SStefano Zampini ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 60877eacf09SStefano Zampini ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr); 60977eacf09SStefano Zampini isbf = (PetscBool)(isbf && v != d); 61077eacf09SStefano Zampini } 61177eacf09SStefano Zampini } 61277eacf09SStefano Zampini if (isbf) { 61377eacf09SStefano Zampini ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr); 61477eacf09SStefano Zampini } 61577eacf09SStefano Zampini } 61677eacf09SStefano Zampini /* count boundary faces */ 61777eacf09SStefano Zampini for (p=fStart,bf=0;p<fEnd;p++) { 61877eacf09SStefano Zampini if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 61977eacf09SStefano Zampini const PetscInt *support; 62077eacf09SStefano Zampini PetscInt supportSize,c; 6218135c375SStefano Zampini 6228135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 6238135c375SStefano Zampini ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 62477eacf09SStefano Zampini if (ovl) { 62577eacf09SStefano Zampini for (c=0;c<supportSize;c++) { 62677eacf09SStefano Zampini if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) { 62777eacf09SStefano Zampini bf++; 6288135c375SStefano Zampini } 62977eacf09SStefano Zampini } 63077eacf09SStefano Zampini } else bf += supportSize; 6318135c375SStefano Zampini } 6328135c375SStefano Zampini } 6338135c375SStefano Zampini /* TODO: is this the label we want to use for marking boundary facets? 6348135c375SStefano Zampini We should decide to have a single marker for these stuff 6358135c375SStefano Zampini Proposal: DMSetBoundaryLabel? 6368135c375SStefano Zampini */ 6378135c375SStefano Zampini ierr = DMGetLabel(dm,"marker",&label);CHKERRQ(ierr); 6388135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr); 6398135c375SStefano Zampini for (p=fStart;p<fEnd;p++) { 64077eacf09SStefano Zampini if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 6418135c375SStefano Zampini const PetscInt *support; 64277eacf09SStefano Zampini PetscInt supportSize,c,nc = 0; 6438135c375SStefano Zampini 6448135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 6458135c375SStefano Zampini ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 64677eacf09SStefano Zampini if (ovl) { 64777eacf09SStefano Zampini for (c=0;c<supportSize;c++) { 64877eacf09SStefano Zampini if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) { 64977eacf09SStefano Zampini fcells[nc++] = support[c]; 6508135c375SStefano Zampini } 65177eacf09SStefano Zampini } 65277eacf09SStefano Zampini } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c]; 65377eacf09SStefano Zampini for (c=0;c<nc;c++) { 65477eacf09SStefano Zampini const PetscInt *cone; 65577eacf09SStefano Zampini int vids[8]; 65677eacf09SStefano Zampini PetscInt i,cell,cl,nv,cid = -1,mid = -1; 6578135c375SStefano Zampini 65877eacf09SStefano Zampini cell = fcells[c]; 65977eacf09SStefano Zampini ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 66077eacf09SStefano Zampini for (cl=0;cl<fpc;cl++) 66177eacf09SStefano Zampini if (cone[cl] == p) 6628135c375SStefano Zampini break; 6638135c375SStefano Zampini 66477eacf09SStefano Zampini /* face material id and type */ 6658135c375SStefano Zampini ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,p,&mid,&cid);CHKERRQ(ierr); 6668135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 66777eacf09SStefano Zampini /* vertex ids */ 66877eacf09SStefano Zampini ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 6698135c375SStefano Zampini for (i=0;i<vpf;i++) { 67077eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer," %D",vids[faces[cl*vpf+i]]);CHKERRQ(ierr); 6718135c375SStefano Zampini } 6728135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 67377eacf09SStefano Zampini } 67477eacf09SStefano Zampini bf = bf-nc; 6758135c375SStefano Zampini } 6768135c375SStefano Zampini } 67777eacf09SStefano Zampini if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf); 67877eacf09SStefano Zampini ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr); 67977eacf09SStefano Zampini ierr = PetscFree(fcells);CHKERRQ(ierr); 6808135c375SStefano Zampini } 6818135c375SStefano Zampini 6828135c375SStefano Zampini /* mark owned vertices */ 6838135c375SStefano Zampini if (ovl) { 6848135c375SStefano Zampini ierr = PetscBTMemzero(vEnd-vStart,pown);CHKERRQ(ierr); 6858135c375SStefano Zampini for (p=cStart;p<cEnd;p++) { 6868135c375SStefano Zampini PetscInt i,closureSize,*closure = NULL; 6878135c375SStefano Zampini 6888135c375SStefano Zampini if (gNum[p-cStart] < 0) continue; 6898135c375SStefano Zampini ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 6908135c375SStefano Zampini for (i=0;i<closureSize;i++) { 6918135c375SStefano Zampini const PetscInt pp = closure[2*i]; 6928135c375SStefano Zampini 6938135c375SStefano Zampini if (pp >= vStart && pp < vEnd) { 6948135c375SStefano Zampini ierr = PetscBTSet(pown,pp-vStart);CHKERRQ(ierr); 6958135c375SStefano Zampini } 6968135c375SStefano Zampini } 6978135c375SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 6988135c375SStefano Zampini } 6998135c375SStefano Zampini } 7008135c375SStefano Zampini if (globalNum) { 7018135c375SStefano Zampini ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 7028135c375SStefano Zampini } 7038135c375SStefano Zampini 7048135c375SStefano Zampini /* vertex_parents (Non-conforming meshes) */ 7058135c375SStefano Zampini parentSection = NULL; 7068135c375SStefano Zampini if (enable_ncmesh) { 7078135c375SStefano Zampini ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 7088135c375SStefano Zampini } 7098135c375SStefano Zampini if (parentSection) { 7108135c375SStefano Zampini PetscInt vp,gvp; 7118135c375SStefano Zampini 7128135c375SStefano Zampini for (vp=0,p=vStart;p<vEnd;p++) { 7138135c375SStefano Zampini DMLabel dlabel; 7148135c375SStefano Zampini PetscInt parent,depth; 7158135c375SStefano Zampini 7168135c375SStefano Zampini if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue; 7178135c375SStefano Zampini ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 7188135c375SStefano Zampini ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr); 7198135c375SStefano Zampini ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 7208135c375SStefano Zampini if (parent != p) vp++; 7218135c375SStefano Zampini } 7228135c375SStefano Zampini ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 7238135c375SStefano Zampini if (gvp) { 7248135c375SStefano Zampini PetscInt maxsupp; 7258135c375SStefano Zampini PetscBool *skip = NULL; 7268135c375SStefano Zampini 7278135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr); 7288135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr); 7298135c375SStefano Zampini ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr); 7308135c375SStefano Zampini ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr); 7318135c375SStefano Zampini for (p=vStart;p<vEnd;p++) { 7328135c375SStefano Zampini DMLabel dlabel; 7338135c375SStefano Zampini PetscInt parent; 7348135c375SStefano Zampini 7358135c375SStefano Zampini if (PetscUnlikely(ovl && !PetscBTLookup(pown,p-vStart))) continue; 7368135c375SStefano Zampini ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 7378135c375SStefano Zampini ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 7388135c375SStefano Zampini if (parent != p) { 739d42aff11SStefano Zampini int vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */ 74011a4995dSStefano Zampini PetscInt i,nv,size,n,numChildren,depth = -1; 7418135c375SStefano Zampini const PetscInt *children; 7428135c375SStefano Zampini ierr = DMPlexGetConeSize(dm,parent,&size);CHKERRQ(ierr); 7438135c375SStefano Zampini switch (size) { 7448135c375SStefano Zampini case 2: /* edge */ 7458135c375SStefano Zampini nv = 0; 746cc0d3ed7SStefano Zampini ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 74777eacf09SStefano Zampini ierr = DMPlexInvertCell(dim,nv,vids);CHKERRQ(ierr); 7488135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr); 7498135c375SStefano Zampini for (i=0;i<nv;i++) { 75011a4995dSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]);CHKERRQ(ierr); 7518135c375SStefano Zampini } 7528135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 7538135c375SStefano Zampini vp--; 7548135c375SStefano Zampini break; 7558135c375SStefano Zampini case 4: /* face */ 7568135c375SStefano Zampini ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr); 7578135c375SStefano Zampini for (n=0;n<numChildren;n++) { 7588135c375SStefano Zampini ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr); 7598135c375SStefano Zampini if (!depth) { 7608135c375SStefano Zampini const PetscInt *hvsupp,*hesupp,*cone; 7618135c375SStefano Zampini PetscInt hvsuppSize,hesuppSize,coneSize; 762451a39c7SStefano Zampini PetscInt hv = children[n],he = -1,f; 7638135c375SStefano Zampini 764451a39c7SStefano Zampini ierr = PetscMemzero(skip,maxsupp*sizeof(PetscBool));CHKERRQ(ierr); 7658135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr); 7668135c375SStefano Zampini ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr); 7678135c375SStefano Zampini for (i=0;i<hvsuppSize;i++) { 7688135c375SStefano Zampini PetscInt ep; 7698135c375SStefano Zampini ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr); 7708135c375SStefano Zampini if (ep != hvsupp[i]) { 7718135c375SStefano Zampini he = hvsupp[i]; 7728135c375SStefano Zampini } else { 7738135c375SStefano Zampini skip[i] = PETSC_TRUE; 7748135c375SStefano Zampini } 7758135c375SStefano Zampini } 776451a39c7SStefano Zampini if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize); 7778135c375SStefano Zampini ierr = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr); 77811a4995dSStefano Zampini vids[0] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 7798135c375SStefano Zampini ierr = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr); 7808135c375SStefano Zampini ierr = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr); 7818135c375SStefano Zampini for (f=0;f<hesuppSize;f++) { 7828135c375SStefano Zampini PetscInt j; 7838135c375SStefano Zampini 7848135c375SStefano Zampini ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr); 7858135c375SStefano Zampini ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr); 7868135c375SStefano Zampini for (j=0;j<coneSize;j++) { 7878135c375SStefano Zampini PetscInt k; 7888135c375SStefano Zampini for (k=0;k<hvsuppSize;k++) { 7898135c375SStefano Zampini if (hvsupp[k] == cone[j]) { 7908135c375SStefano Zampini skip[k] = PETSC_TRUE; 7918135c375SStefano Zampini break; 7928135c375SStefano Zampini } 7938135c375SStefano Zampini } 7948135c375SStefano Zampini } 7958135c375SStefano Zampini } 7968135c375SStefano Zampini for (i=0;i<hvsuppSize;i++) { 7978135c375SStefano Zampini if (!skip[i]) { 7988135c375SStefano Zampini ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr); 79911a4995dSStefano Zampini vids[1] = (int)((cone[0] == hv) ? cone[1] : cone[0]); 8008135c375SStefano Zampini } 8018135c375SStefano Zampini } 8028135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr); 8038135c375SStefano Zampini for (i=0;i<2;i++) { 80411a4995dSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer," %d",vids[i]-vStart);CHKERRQ(ierr); 8058135c375SStefano Zampini } 8068135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 8078135c375SStefano Zampini vp--; 8088135c375SStefano Zampini } 8098135c375SStefano Zampini } 8108135c375SStefano Zampini break; 8118135c375SStefano Zampini default: 8128135c375SStefano Zampini SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %d",size); 8138135c375SStefano Zampini } 8148135c375SStefano Zampini } 8158135c375SStefano Zampini } 8168135c375SStefano Zampini ierr = PetscFree(skip);CHKERRQ(ierr); 8178135c375SStefano Zampini } 8188135c375SStefano Zampini if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp); 8198135c375SStefano Zampini } 8208135c375SStefano Zampini ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 8218135c375SStefano Zampini 8228135c375SStefano Zampini /* vertices */ 82377eacf09SStefano Zampini if (hovec) { /* higher-order meshes */ 82477eacf09SStefano Zampini const char *fec; 82577eacf09SStefano Zampini PetscInt i,n; 82677eacf09SStefano Zampini 82777eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 82877eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr); 82977eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr); 83077eacf09SStefano Zampini ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr); 83177eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr); 83277eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementCollection: %s\n",fec);CHKERRQ(ierr); 83377eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr); 83477eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/ 83577eacf09SStefano Zampini ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr); 83677eacf09SStefano Zampini ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr); 83777eacf09SStefano Zampini if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim); 83877eacf09SStefano Zampini for (i=0;i<n/sdim;i++) { 83977eacf09SStefano Zampini PetscInt s; 84077eacf09SStefano Zampini for (s=0;s<sdim;s++) { 84177eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr); 84277eacf09SStefano Zampini } 84377eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 84477eacf09SStefano Zampini } 84577eacf09SStefano Zampini ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr); 84677eacf09SStefano Zampini } else { 847cc0d3ed7SStefano Zampini ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr); 8488135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 849cc0d3ed7SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr); 8508135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 8518135c375SStefano Zampini ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr); 852cc0d3ed7SStefano Zampini for (p=0;p<nvert/sdim;p++) { 853cc0d3ed7SStefano Zampini PetscInt s; 854cc0d3ed7SStefano Zampini for (s=0;s<sdim;s++) { 85577eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[p*sdim+s]));CHKERRQ(ierr); 8568135c375SStefano Zampini } 857cc0d3ed7SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 8588135c375SStefano Zampini } 8598135c375SStefano Zampini ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr); 86077eacf09SStefano Zampini } 8618135c375SStefano Zampini PetscFunctionReturn(0); 8628135c375SStefano Zampini } 8638135c375SStefano Zampini 8648135c375SStefano Zampini /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump */ 8658135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer) 8668135c375SStefano Zampini { 8678135c375SStefano Zampini PetscErrorCode ierr; 8688135c375SStefano Zampini PetscBool isglvis,isascii; 8698135c375SStefano Zampini 8708135c375SStefano Zampini PetscFunctionBegin; 8718135c375SStefano Zampini PetscValidHeaderSpecific(dm,DM_CLASSID,1); 8728135c375SStefano Zampini PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 8738135c375SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr); 8748135c375SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 8758135c375SStefano Zampini if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII"); 8768135c375SStefano Zampini if (isglvis) { 8778135c375SStefano Zampini PetscViewer view; 8788135c375SStefano Zampini PetscViewerGLVisType type; 8798135c375SStefano Zampini 8808135c375SStefano Zampini ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr); 8818135c375SStefano Zampini ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr); 8828135c375SStefano Zampini if (view) { /* in the socket case, it may happen that the connection failed */ 8838135c375SStefano Zampini if (type == PETSC_VIEWER_GLVIS_SOCKET) { 8848135c375SStefano Zampini PetscMPIInt size,rank; 8858135c375SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 8868135c375SStefano Zampini ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 8878135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr); 8888135c375SStefano Zampini } 8898135c375SStefano Zampini ierr = DMPlexView_GLVis_ASCII(dm,view);CHKERRQ(ierr); 8908135c375SStefano Zampini ierr = PetscViewerFlush(view);CHKERRQ(ierr); 8918135c375SStefano Zampini if (type == PETSC_VIEWER_GLVIS_SOCKET) { 8928135c375SStefano Zampini PetscInt dim; 8938135c375SStefano Zampini const char* name; 8948135c375SStefano Zampini 8958135c375SStefano Zampini ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr); 8968135c375SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 8978135c375SStefano Zampini ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr); 8988135c375SStefano Zampini ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr); 8998135c375SStefano Zampini } 9008135c375SStefano Zampini } 9018135c375SStefano Zampini } else { 9028135c375SStefano Zampini ierr = DMPlexView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr); 9038135c375SStefano Zampini } 9048135c375SStefano Zampini PetscFunctionReturn(0); 9058135c375SStefano Zampini } 906