#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; PetscErrorCode ierr; PetscFunctionBegin; for (i=0;inf;i++) { ierr = VecScatterDestroy(&ctx->scctx[i]);CHKERRQ(ierr); } ierr = PetscFree(ctx->scctx);CHKERRQ(ierr); ierr = PetscFree(vctx);CHKERRQ(ierr); PetscFunctionReturn(0); } static PetscErrorCode DMPlexSampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXfield[], void *vctx) { GLVisViewerCtx *ctx = (GLVisViewerCtx*)vctx; PetscInt f; PetscErrorCode ierr; PetscFunctionBegin; for (f=0;fscctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(ctx->scctx[f],(Vec)oX,(Vec)oXfield[f],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); } 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; PetscErrorCode ierr; PetscFunctionBegin; ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr); if (!globalNum) { ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr); ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr); ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr); } ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); for (c = cStart, totc = 0; c < cEnd; c++) { if (gNum[c-cStart] >= 0) { PetscInt i,numPoints,*points = NULL; totc++; ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); for (i=0;i= vStart) && (points[i] < vEnd)) { ierr = PetscBTSet(vown,points[i]-vStart);CHKERRQ(ierr); } } ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); } } for (f=0,Nv=0;fscctx);CHKERRQ(ierr); ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); if (ds) { for (f=0;f 0) { /* no support for high-order viz, still have to figure out the numbering */ ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim);CHKERRQ(ierr); } else { if (!continuous && order) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported for order %D",order); H1 = PETSC_FALSE; ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P%D",dim,order);CHKERRQ(ierr); } ierr = PetscStrallocpy(name,&fieldname[ctx->nf]);CHKERRQ(ierr); bs[ctx->nf] = Nc; dims[ctx->nf] = dim; if (H1) { nlocal[ctx->nf] = Nc * Nv; ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr); for (i=0,cum=0;inf] = Nc * totc; ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr); for (i=0,cum=0;iscctx[ctx->nf]);CHKERRQ(ierr); ierr = VecDestroy(&xfield);CHKERRQ(ierr); ierr = ISDestroy(&isfield);CHKERRQ(ierr); ctx->nf++; } else if (id == PETSCFV_CLASSID) { PetscInt c; ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr); ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P0",dim);CHKERRQ(ierr); for (c = 0; c < Nc; c++) { char comp[256]; ierr = PetscSNPrintf(comp,256,"%s-Comp%D",name,c);CHKERRQ(ierr); ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr); bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */ nlocal[ctx->nf] = totc; dims[ctx->nf] = dim; ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr); for (i=0,cum=0;iscctx[ctx->nf]);CHKERRQ(ierr); ierr = VecDestroy(&xfield);CHKERRQ(ierr); ierr = ISDestroy(&isfield);CHKERRQ(ierr); ctx->nf++; } } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f); } else SETERRQ1(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"); ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); ierr = VecDestroy(&xlocal);CHKERRQ(ierr); ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); /* create work vectors */ for (f=0;fnf;f++) { ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]);CHKERRQ(ierr); ierr = VecSetBlockSize(Ufield[f],bs[f]);CHKERRQ(ierr); ierr = VecSetDM(Ufield[f],dm);CHKERRQ(ierr); } /* customize the viewer */ ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr); for (f=0;fnf;f++) { ierr = PetscFree(fieldname[f]);CHKERRQ(ierr); ierr = PetscFree(fec_type[f]);CHKERRQ(ierr); ierr = VecDestroy(&Ufield[f]);CHKERRQ(ierr); } ierr = PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield);CHKERRQ(ierr); 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; PetscErrorCode ierr; PetscFunctionBegin; ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); ierr = DMLabelGetValue(dlabel,p,&pdepth);CHKERRQ(ierr); ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr); ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); if (label) { ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr); *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 */ if (dim < 0 || dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim); if (csize > 8) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p); if (depth != 1) SETERRQ2(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 { if (csize > 6) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p); if (pdepth < 0 || pdepth > 3) SETERRQ2(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; PetscErrorCode ierr; PetscFunctionBegin; ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); sdim = dim; if (csec) { PetscInt sStart,sEnd; ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); ierr = PetscSectionGetChart(csec,&sStart,&sEnd);CHKERRQ(ierr); ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr); off = off/sdim; if (p >= sStart && p < sEnd) { ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); } } if (!dof) { ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); for (i=0,q=0;i= vStart) && (points[i] < vEnd)) vids[q++] = points[i]-vStart+off; ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); } else { ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr); ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); for (q=0;q 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); /* get container: determines if a process visualizes is portion of the data or not */ ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr); if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); { PetscViewerGLVisInfo glvis_info; ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr); 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 DMPlex does not currently support HO meshes, so there's no API for this */ ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr); ierr = PetscObjectReference((PetscObject)hovec);CHKERRQ(ierr); if (!hovec) { DM cdm; PetscFE disc; PetscClassId classid; ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); ierr = DMGetField(cdm,0,NULL,(PetscObject*)&disc);CHKERRQ(ierr); ierr = PetscObjectGetClassId((PetscObject)disc,&classid);CHKERRQ(ierr); if (classid == PETSCFE_CLASSID) { DM hocdm; PetscFE hodisc; Vec vec; Mat mat; char name[32],fec_type[64]; ierr = GLVisCreateFE(disc,name,&hodisc);CHKERRQ(ierr); ierr = DMClone(cdm,&hocdm);CHKERRQ(ierr); ierr = DMSetField(hocdm,0,NULL,(PetscObject)hodisc);CHKERRQ(ierr); ierr = PetscFEDestroy(&hodisc);CHKERRQ(ierr); ierr = DMCreateDS(hocdm);CHKERRQ(ierr); ierr = DMGetCoordinates(dm,&vec);CHKERRQ(ierr); ierr = DMCreateGlobalVector(hocdm,&hovec);CHKERRQ(ierr); ierr = DMCreateInterpolation(cdm,hocdm,&mat,NULL);CHKERRQ(ierr); ierr = MatInterpolate(mat,vec,hovec);CHKERRQ(ierr); ierr = MatDestroy(&mat);CHKERRQ(ierr); ierr = DMDestroy(&hocdm);CHKERRQ(ierr); ierr = PetscSNPrintf(fec_type,sizeof(fec_type),"FiniteElementCollection: %s", name);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject)hovec,fec_type);CHKERRQ(ierr); } } ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); ierr = DMPlexGetGhostCellStratum(dm,&p,NULL);CHKERRQ(ierr); if (p >= 0) cEnd = p; ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr); ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr); if (periodic && !localized && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized"); ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr); if (!coordinates && !hovec) SETERRQ(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");CHKERRQ(ierr); ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr); 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); 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); ierr = PetscOptionsBool("-viewer_glvis_dm_plex_overlap","Include overlap region in local meshes",NULL,view_ovl,&view_ovl,NULL);CHKERRQ(ierr); ierr = PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark);CHKERRQ(ierr); ierr = PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark);CHKERRQ(ierr); ierr = PetscOptionsEnd();CHKERRQ(ierr); if (enable_bmark) enable_boundary = PETSC_TRUE; ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRMPI(ierr); if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel"); ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0"); if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(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 */ if (depth != 1) SETERRQ1(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; ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr); if (!globalNum) { if (view_ovl) { ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),cEnd-cStart,0,1,&globalNum);CHKERRQ(ierr); } else { ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr); } ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr); ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr); } ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); for (p=cStart; p= vStart && cellClosure[2*v] < vEnd) { vidxs = cellClosure + 2*v; break; } if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices"); for (v=0;vPETSC_MACHINE_EPSILON) { ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr); } } } ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); } } if (dim > 1) { PetscInt eEnd,eStart; ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr); 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; ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); isbf = PETSC_TRUE; for (i=0;i= vStart && pp < vEnd) { ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr); } } ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); } } /* vertex_parents (Non-conforming meshes) */ parentSection = NULL; if (enable_ncmesh) { ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); } if (parentSection) { PetscInt vp,gvp; for (vp=0,p=vStart;p