#include #include #include #include #include #include #include typedef struct { PetscInt nf; VecScatter *scctx; } GLVisViewerCtx; static PetscErrorCode DestroyGLVisViewerCtx_Private(void *vctx) { GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; PetscInt i; PetscFunctionBegin; for (i=0;inf;i++) { PetscCall(VecScatterDestroy(&ctx->scctx[i])); } PetscCall(PetscFree(ctx->scctx)); PetscCall(PetscFree(vctx)); PetscFunctionReturn(0); } static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) { GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; PetscInt f; PetscFunctionBegin; for (f=0;fscctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD)); PetscCall(VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD)); } PetscFunctionReturn(0); } /* for FEM, it works for H1 fields only and extracts dofs at cell vertices, discarding any other dof */ PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject odm, PetscViewer viewer) { DM dm = (DM)odm; Vec xlocal,xfield,*Ufield; PetscDS ds; IS globalNum,isfield; PetscBT vown; char **fieldname = NULL,**fec_type = NULL; const PetscInt *gNum; PetscInt *nlocal,*bs,*idxs,*dims; PetscInt f,maxfields,nfields,c,totc,totdofs,Nv,cum,i; PetscInt dim,cStart,cEnd,vStart,vEnd; GLVisViewerCtx *ctx; PetscSection s; PetscFunctionBegin; PetscCall(DMGetDimension(dm,&dim)); PetscCall(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd)); PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd)); PetscCall(PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum)); if (!globalNum) { PetscCall(DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum)); PetscCall(PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum)); PetscCall(PetscObjectDereference((PetscObject)globalNum)); } PetscCall(ISGetIndices(globalNum,&gNum)); PetscCall(PetscBTCreate(vEnd-vStart,&vown)); for (c = cStart, totc = 0; c < cEnd; c++) { if (gNum[c-cStart] >= 0) { PetscInt i,numPoints,*points = NULL; totc++; PetscCall(DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points)); for (i=0;i= vStart) && (points[i] < vEnd)) { PetscCall(PetscBTSet(vown,points[i]-vStart)); } } PetscCall(DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points)); } } for (f=0,Nv=0;fscctx)); PetscCall(DMGetDS(dm,&ds)); if (ds) { for (f=0;f 0) { /* no support for high-order viz, still have to figure out the numbering */ PetscCall(PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim)); } else { PetscCheckFalse(!continuous && order,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported for order %D",order); H1 = PETSC_FALSE; PetscCall(PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P%D",dim,order)); } PetscCall(PetscStrallocpy(name,&fieldname[ctx->nf])); bs[ctx->nf] = Nc; dims[ctx->nf] = dim; if (H1) { nlocal[ctx->nf] = Nc * Nv; PetscCall(PetscStrallocpy(fec,&fec_type[ctx->nf])); PetscCall(VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield)); for (i=0,cum=0;inf] = Nc * totc; PetscCall(PetscStrallocpy(fec,&fec_type[ctx->nf])); PetscCall(VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield)); for (i=0,cum=0;iscctx[ctx->nf])); PetscCall(VecDestroy(&xfield)); PetscCall(ISDestroy(&isfield)); ctx->nf++; } else if (id == PETSCFV_CLASSID) { PetscInt c; PetscCall(PetscFVGetNumComponents((PetscFV)disc,&Nc)); PetscCall(PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P0",dim)); for (c = 0; c < Nc; c++) { char comp[256]; PetscCall(PetscSNPrintf(comp,256,"%s-Comp%D",name,c)); PetscCall(PetscStrallocpy(comp,&fieldname[ctx->nf])); bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */ nlocal[ctx->nf] = totc; dims[ctx->nf] = dim; PetscCall(PetscStrallocpy(fec,&fec_type[ctx->nf])); PetscCall(VecCreateSeq(PETSC_COMM_SELF,totc,&xfield)); for (i=0,cum=0;iscctx[ctx->nf])); PetscCall(VecDestroy(&xfield)); PetscCall(ISDestroy(&isfield)); ctx->nf++; } } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f); } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f); } } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM"); PetscCall(PetscBTDestroy(&vown)); PetscCall(VecDestroy(&xlocal)); PetscCall(ISRestoreIndices(globalNum,&gNum)); /* create work vectors */ for (f=0;fnf;f++) { PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f])); PetscCall(PetscObjectSetName((PetscObject)Ufield[f],fieldname[f])); PetscCall(VecSetBlockSize(Ufield[f],bs[f])); PetscCall(VecSetDM(Ufield[f],dm)); } /* customize the viewer */ PetscCall(PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private)); for (f=0;fnf;f++) { PetscCall(PetscFree(fieldname[f])); PetscCall(PetscFree(fec_type[f])); PetscCall(VecDestroy(&Ufield[f])); } PetscCall(PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield)); PetscFunctionReturn(0); } typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_PRISM,MFEM_UNDEF} MFEM_cid; MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_UNDEF}, {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_UNDEF}, {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF,MFEM_UNDEF}, {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_PRISM,MFEM_CUBE } }; MFEM_cid mfem_table_cid_unint[4][9] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_CUBE } }; static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt minl, PetscInt p, PetscInt *mid, PetscInt *cid) { DMLabel dlabel; PetscInt depth,csize,pdepth,dim; PetscFunctionBegin; PetscCall(DMPlexGetDepthLabel(dm,&dlabel)); PetscCall(DMLabelGetValue(dlabel,p,&pdepth)); PetscCall(DMPlexGetConeSize(dm,p,&csize)); PetscCall(DMPlexGetDepth(dm,&depth)); PetscCall(DMGetDimension(dm,&dim)); if (label) { PetscCall(DMLabelGetValue(label,p,mid)); *mid = *mid - minl + 1; /* MFEM does not like negative markers */ } else *mid = 1; if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ PetscCheckFalse(dim < 0 || dim > 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim); PetscCheckFalse(csize > 8,PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p); PetscCheckFalse(depth != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Found depth %D for point %D. You should interpolate the mesh first",depth,p); *cid = mfem_table_cid_unint[dim][csize]; } else { PetscCheckFalse(csize > 6,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p); PetscCheckFalse(pdepth < 0 || pdepth > 3,PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p); *cid = mfem_table_cid[pdepth][csize]; } PetscFunctionReturn(0); } static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, PetscInt vids[]) { PetscInt dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL; PetscFunctionBegin; PetscCall(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd)); PetscCall(DMGetDimension(dm,&dim)); sdim = dim; if (csec) { PetscInt sStart,sEnd; PetscCall(DMGetCoordinateDim(dm,&sdim)); PetscCall(PetscSectionGetChart(csec,&sStart,&sEnd)); PetscCall(PetscSectionGetOffset(csec,vStart,&off)); off = off/sdim; if (p >= sStart && p < sEnd) { PetscCall(PetscSectionGetDof(csec,p,&dof)); } } if (!dof) { PetscCall(DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points)); for (i=0,q=0;i= vStart) && (points[i] < vEnd)) vids[q++] = points[i]-vStart+off; PetscCall(DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points)); } else { PetscCall(PetscSectionGetOffset(csec,p,&off)); PetscCall(PetscSectionGetDof(csec,p,&dof)); for (q=0;q 1,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); PetscCall(DMGetDimension(dm,&dim)); PetscCall(DMPlexGetDepth(dm,&depth)); /* get container: determines if a process visualizes is portion of the data or not */ PetscCall(PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container)); PetscCheck(glvis_container,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); { PetscViewerGLVisInfo glvis_info; PetscCall(PetscContainerGetPointer(glvis_container,(void**)&glvis_info)); enabled = glvis_info->enabled; fmt = glvis_info->fmt; } /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh */ PetscCall(PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec)); PetscCall(PetscObjectReference((PetscObject)hovec)); if (!hovec) { DM cdm; PetscFE disc; PetscClassId classid; PetscCall(DMGetCoordinateDM(dm,&cdm)); PetscCall(DMGetField(cdm,0,NULL,(PetscObject*)&disc)); PetscCall(PetscObjectGetClassId((PetscObject)disc,&classid)); if (classid == PETSCFE_CLASSID) { DM hocdm; PetscFE hodisc; Vec vec; Mat mat; char name[32],fec_type[64]; IS perm = NULL; PetscCall(GLVisCreateFE(disc,name,&hodisc,&perm)); PetscCall(DMClone(cdm,&hocdm)); PetscCall(DMSetField(hocdm,0,NULL,(PetscObject)hodisc)); PetscCall(PetscFEDestroy(&hodisc)); PetscCall(DMCreateDS(hocdm)); PetscCall(DMGetCoordinates(dm,&vec)); PetscCall(DMCreateGlobalVector(hocdm,&hovec)); PetscCall(DMCreateInterpolation(cdm,hocdm,&mat,NULL)); PetscCall(MatInterpolate(mat,vec,hovec)); PetscCall(MatDestroy(&mat)); PetscCall(DMGetLocalSection(hocdm,&hoSection)); PetscCall(PetscSectionSetClosurePermutation(hoSection, (PetscObject)hocdm, depth, perm)); PetscCall(ISDestroy(&perm)); PetscCall(DMDestroy(&hocdm)); PetscCall(PetscSNPrintf(fec_type,sizeof(fec_type),"FiniteElementCollection: %s", name)); PetscCall(PetscObjectSetName((PetscObject)hovec,fec_type)); } } PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd)); PetscCall(DMPlexGetGhostCellStratum(dm,&p,NULL)); if (p >= 0) cEnd = p; PetscCall(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd)); PetscCall(DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL)); PetscCall(DMGetCoordinatesLocalized(dm,&localized)); PetscCheckFalse(periodic && !localized && !hovec,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized"); PetscCall(DMGetCoordinateSection(dm,&coordSection)); PetscCall(DMGetCoordinateDim(dm,&sdim)); PetscCall(DMGetCoordinatesLocal(dm,&coordinates)); PetscCheckFalse(!coordinates && !hovec,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector"); /* a couple of sections of the mesh specification are disabled - boundary: the boundary is not needed for proper mesh visualization unless we want to visualize boundary attributes or we have high-order coordinates in 3D (topologically) - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel) */ enable_boundary = PETSC_FALSE; enable_ncmesh = PETSC_FALSE; enable_mfem = PETSC_FALSE; enable_emark = PETSC_FALSE; enable_bmark = PETSC_FALSE; /* I'm tired of problems with negative values in the markers, disable them */ ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");PetscCall(ierr); PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL)); PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_ncmesh","Enable vertex_parents section in mesh representation (allows derefinement)",NULL,enable_ncmesh,&enable_ncmesh,NULL)); PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_enable_mfem","Dump a mesh that can be used with MFEM's FiniteElementSpaces",NULL,enable_mfem,&enable_mfem,NULL)); PetscCall(PetscOptionsBool("-viewer_glvis_dm_plex_overlap","Include overlap region in local meshes",NULL,view_ovl,&view_ovl,NULL)); PetscCall(PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark)); PetscCall(PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark)); ierr = PetscOptionsEnd();PetscCall(ierr); if (enable_bmark) enable_boundary = PETSC_TRUE; PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size)); PetscCheckFalse(enable_ncmesh && size > 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel"); PetscCheckFalse(enable_boundary && depth >= 0 && dim != depth,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0"); PetscCheckFalse(enable_ncmesh && depth >= 0 && dim != depth,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0"); if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ PetscCheckFalse(depth != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth); cellvertex = PETSC_TRUE; } /* Identify possible cells in the overlap */ novl = 0; pown = NULL; if (size > 1) { IS globalNum = NULL; const PetscInt *gNum; PetscBool ovl = PETSC_FALSE; PetscCall(PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum)); if (!globalNum) { if (view_ovl) { PetscCall(ISCreateStride(PetscObjectComm((PetscObject)dm),cEnd-cStart,0,1,&globalNum)); } else { PetscCall(DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum)); } PetscCall(PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum)); PetscCall(PetscObjectDereference((PetscObject)globalNum)); } PetscCall(ISGetIndices(globalNum,&gNum)); for (p=cStart; p= vStart && cellClosure[2*v] < vEnd) { vidxs = cellClosure + 2*v; break; } PetscCheck(vidxs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices"); for (v=0;vPETSC_MACHINE_EPSILON) { PetscCall(DMLabelSetValue(perLabel,vidxs[2*v],2)); } } } PetscCall(DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure)); PetscCall(DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals)); } } if (dim > 1) { PetscInt eEnd,eStart; PetscCall(DMPlexGetDepthStratum(dm,1,&eStart,&eEnd)); for (p=eStart;p 2) { for (p=fStart;p 1 && has_owned && has_ghost)); } else { isbf = (PetscBool)(supportSize == 1); } if (!isbf && perLabel) { const PetscInt *cone; PetscInt coneSize,i; PetscCall(DMPlexGetCone(dm,p,&cone)); PetscCall(DMPlexGetConeSize(dm,p,&coneSize)); isbf = PETSC_TRUE; for (i=0;i= vStart && pp < vEnd) { PetscCall(PetscBTSet(vown,pp-vStart)); } } PetscCall(DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure)); } } if (parentSection) { PetscInt vp,gvp; for (vp=0,p=vStart;p