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,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 = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr); 65 if (!globalNum) { 66 ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr); 67 ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr); 68 ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr); 69 } 70 ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 71 ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); 72 for (c = cStart, totc = 0; c < cEnd; c++) { 73 if (gNum[c-cStart] >= 0) { 74 PetscInt i,numPoints,*points = NULL; 75 76 totc++; 77 ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 78 for (i=0;i<numPoints*2;i+= 2) { 79 if ((points[i] >= vStart) && (points[i] < vEnd)) { 80 ierr = PetscBTSet(vown,points[i]-vStart);CHKERRQ(ierr); 81 } 82 } 83 ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 84 } 85 } 86 for (f=0,Nv=0;f<vEnd-vStart;f++) if (PetscLikely(PetscBTLookup(vown,f))) Nv++; 87 88 ierr = DMCreateLocalVector(dm,&xlocal);CHKERRQ(ierr); 89 ierr = VecGetLocalSize(xlocal,&totdofs);CHKERRQ(ierr); 90 ierr = DMGetLocalSection(dm,&s);CHKERRQ(ierr); 91 ierr = PetscSectionGetNumFields(s,&nfields);CHKERRQ(ierr); 92 for (f=0,maxfields=0;f<nfields;f++) { 93 PetscInt bs; 94 95 ierr = PetscSectionGetFieldComponents(s,f,&bs);CHKERRQ(ierr); 96 maxfields += bs; 97 } 98 ierr = PetscCalloc7(maxfields,&fieldname,maxfields,&nlocal,maxfields,&bs,maxfields,&dims,maxfields,&fec_type,totdofs,&idxs,maxfields,&Ufield);CHKERRQ(ierr); 99 ierr = PetscNew(&ctx);CHKERRQ(ierr); 100 ierr = PetscCalloc1(maxfields,&ctx->scctx);CHKERRQ(ierr); 101 ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 102 if (ds) { 103 for (f=0;f<nfields;f++) { 104 const char* fname; 105 char name[256]; 106 PetscObject disc; 107 size_t len; 108 109 ierr = PetscSectionGetFieldName(s,f,&fname);CHKERRQ(ierr); 110 ierr = PetscStrlen(fname,&len);CHKERRQ(ierr); 111 if (len) { 112 ierr = PetscStrcpy(name,fname);CHKERRQ(ierr); 113 } else { 114 ierr = PetscSNPrintf(name,256,"Field%D",f);CHKERRQ(ierr); 115 } 116 ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr); 117 if (disc) { 118 PetscClassId id; 119 PetscInt Nc; 120 char fec[64]; 121 122 ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr); 123 if (id == PETSCFE_CLASSID) { 124 PetscFE fem = (PetscFE)disc; 125 PetscDualSpace sp; 126 PetscDualSpaceType spname; 127 PetscInt order; 128 PetscBool islag,continuous,H1 = PETSC_TRUE; 129 130 ierr = PetscFEGetNumComponents(fem,&Nc);CHKERRQ(ierr); 131 ierr = PetscFEGetDualSpace(fem,&sp);CHKERRQ(ierr); 132 ierr = PetscDualSpaceGetType(sp,&spname);CHKERRQ(ierr); 133 ierr = PetscStrcmp(spname,PETSCDUALSPACELAGRANGE,&islag);CHKERRQ(ierr); 134 if (!islag) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupported dual space"); 135 ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr); 136 ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr); 137 if (continuous && order > 0) { /* no support for high-order viz, still have to figure out the numbering */ 138 ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: H1_%DD_P1",dim);CHKERRQ(ierr); 139 } else { 140 if (!continuous && order) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Discontinuous space visualization currently unsupported for order %D",order); 141 H1 = PETSC_FALSE; 142 ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P%D",dim,order);CHKERRQ(ierr); 143 } 144 ierr = PetscStrallocpy(name,&fieldname[ctx->nf]);CHKERRQ(ierr); 145 bs[ctx->nf] = Nc; 146 dims[ctx->nf] = dim; 147 if (H1) { 148 nlocal[ctx->nf] = Nc * Nv; 149 ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 150 ierr = VecCreateSeq(PETSC_COMM_SELF,Nv*Nc,&xfield);CHKERRQ(ierr); 151 for (i=0,cum=0;i<vEnd-vStart;i++) { 152 PetscInt j,off; 153 154 if (PetscUnlikely(!PetscBTLookup(vown,i))) continue; 155 ierr = PetscSectionGetFieldOffset(s,i+vStart,f,&off);CHKERRQ(ierr); 156 for (j=0;j<Nc;j++) idxs[cum++] = off + j; 157 } 158 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),Nv*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 159 } else { 160 nlocal[ctx->nf] = Nc * totc; 161 ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 162 ierr = VecCreateSeq(PETSC_COMM_SELF,Nc*totc,&xfield);CHKERRQ(ierr); 163 for (i=0,cum=0;i<cEnd-cStart;i++) { 164 PetscInt j,off; 165 166 if (PetscUnlikely(gNum[i] < 0)) continue; 167 ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr); 168 for (j=0;j<Nc;j++) idxs[cum++] = off + j; 169 } 170 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc*Nc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 171 } 172 ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr); 173 ierr = VecDestroy(&xfield);CHKERRQ(ierr); 174 ierr = ISDestroy(&isfield);CHKERRQ(ierr); 175 ctx->nf++; 176 } else if (id == PETSCFV_CLASSID) { 177 PetscInt c; 178 179 ierr = PetscFVGetNumComponents((PetscFV)disc,&Nc);CHKERRQ(ierr); 180 ierr = PetscSNPrintf(fec,64,"FiniteElementCollection: L2_%DD_P0",dim);CHKERRQ(ierr); 181 for (c = 0; c < Nc; c++) { 182 char comp[256]; 183 ierr = PetscSNPrintf(comp,256,"%s-Comp%D",name,c);CHKERRQ(ierr); 184 ierr = PetscStrallocpy(comp,&fieldname[ctx->nf]);CHKERRQ(ierr); 185 bs[ctx->nf] = 1; /* Does PetscFV support components with different block size? */ 186 nlocal[ctx->nf] = totc; 187 dims[ctx->nf] = dim; 188 ierr = PetscStrallocpy(fec,&fec_type[ctx->nf]);CHKERRQ(ierr); 189 ierr = VecCreateSeq(PETSC_COMM_SELF,totc,&xfield);CHKERRQ(ierr); 190 for (i=0,cum=0;i<cEnd-cStart;i++) { 191 PetscInt off; 192 193 if (PetscUnlikely(gNum[i])<0) continue; 194 ierr = PetscSectionGetFieldOffset(s,i+cStart,f,&off);CHKERRQ(ierr); 195 idxs[cum++] = off + c; 196 } 197 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)xlocal),totc,idxs,PETSC_USE_POINTER,&isfield);CHKERRQ(ierr); 198 ierr = VecScatterCreate(xlocal,isfield,xfield,NULL,&ctx->scctx[ctx->nf]);CHKERRQ(ierr); 199 ierr = VecDestroy(&xfield);CHKERRQ(ierr); 200 ierr = ISDestroy(&isfield);CHKERRQ(ierr); 201 ctx->nf++; 202 } 203 } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONG,"Unknown discretization type for field %D",f); 204 } else SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing discretization for field %D",f); 205 } 206 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Needs a DS attached to the DM"); 207 ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); 208 ierr = VecDestroy(&xlocal);CHKERRQ(ierr); 209 ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 210 211 /* create work vectors */ 212 for (f=0;f<ctx->nf;f++) { 213 ierr = VecCreateMPI(PetscObjectComm((PetscObject)dm),nlocal[f],PETSC_DECIDE,&Ufield[f]);CHKERRQ(ierr); 214 ierr = PetscObjectSetName((PetscObject)Ufield[f],fieldname[f]);CHKERRQ(ierr); 215 ierr = VecSetBlockSize(Ufield[f],bs[f]);CHKERRQ(ierr); 216 ierr = VecSetDM(Ufield[f],dm);CHKERRQ(ierr); 217 } 218 219 /* customize the viewer */ 220 ierr = PetscViewerGLVisSetFields(viewer,ctx->nf,(const char**)fec_type,dims,DMPlexSampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DestroyGLVisViewerCtx_Private);CHKERRQ(ierr); 221 for (f=0;f<ctx->nf;f++) { 222 ierr = PetscFree(fieldname[f]);CHKERRQ(ierr); 223 ierr = PetscFree(fec_type[f]);CHKERRQ(ierr); 224 ierr = VecDestroy(&Ufield[f]);CHKERRQ(ierr); 225 } 226 ierr = PetscFree7(fieldname,nlocal,bs,dims,fec_type,idxs,Ufield);CHKERRQ(ierr); 227 PetscFunctionReturn(0); 228 } 229 230 typedef enum {MFEM_POINT=0,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE,MFEM_TETRAHEDRON,MFEM_CUBE,MFEM_PRISM,MFEM_UNDEF} MFEM_cid; 231 232 MFEM_cid mfem_table_cid[4][7] = { {MFEM_POINT,MFEM_UNDEF,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_UNDEF}, 233 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_UNDEF}, 234 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF,MFEM_UNDEF}, 235 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_PRISM,MFEM_CUBE } }; 236 237 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}, 238 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_UNDEF ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, 239 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_TRIANGLE,MFEM_SQUARE ,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_UNDEF}, 240 {MFEM_POINT,MFEM_UNDEF,MFEM_SEGMENT,MFEM_UNDEF ,MFEM_TETRAHEDRON,MFEM_UNDEF,MFEM_PRISM,MFEM_UNDEF,MFEM_CUBE } }; 241 242 static PetscErrorCode DMPlexGetPointMFEMCellID_Internal(DM dm, DMLabel label, PetscInt minl, PetscInt p, PetscInt *mid, PetscInt *cid) 243 { 244 DMLabel dlabel; 245 PetscInt depth,csize,pdepth,dim; 246 PetscErrorCode ierr; 247 248 PetscFunctionBegin; 249 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 250 ierr = DMLabelGetValue(dlabel,p,&pdepth);CHKERRQ(ierr); 251 ierr = DMPlexGetConeSize(dm,p,&csize);CHKERRQ(ierr); 252 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 253 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 254 if (label) { 255 ierr = DMLabelGetValue(label,p,mid);CHKERRQ(ierr); 256 *mid = *mid - minl + 1; /* MFEM does not like negative markers */ 257 } else *mid = 1; 258 if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 259 if (dim < 0 || dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Dimension %D",dim); 260 if (csize > 8) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found cone size %D for point %D",csize,p); 261 if (depth != 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Found depth %D for point %D. You should interpolate the mesh first",depth,p); 262 *cid = mfem_table_cid_unint[dim][csize]; 263 } else { 264 if (csize > 6) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cone size %D for point %D",csize,p); 265 if (pdepth < 0 || pdepth > 3) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Depth %D for point %D",csize,p); 266 *cid = mfem_table_cid[pdepth][csize]; 267 } 268 PetscFunctionReturn(0); 269 } 270 271 static PetscErrorCode DMPlexGetPointMFEMVertexIDs_Internal(DM dm, PetscInt p, PetscSection csec, PetscInt *nv, PetscInt vids[]) 272 { 273 PetscInt dim,sdim,dof = 0,off = 0,i,q,vStart,vEnd,numPoints,*points = NULL; 274 PetscErrorCode ierr; 275 276 PetscFunctionBegin; 277 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 278 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 279 sdim = dim; 280 if (csec) { 281 PetscInt sStart,sEnd; 282 283 ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 284 ierr = PetscSectionGetChart(csec,&sStart,&sEnd);CHKERRQ(ierr); 285 ierr = PetscSectionGetOffset(csec,vStart,&off);CHKERRQ(ierr); 286 off = off/sdim; 287 if (p >= sStart && p < sEnd) { 288 ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 289 } 290 } 291 if (!dof) { 292 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 293 for (i=0,q=0;i<numPoints*2;i+= 2) 294 if ((points[i] >= vStart) && (points[i] < vEnd)) 295 vids[q++] = points[i]-vStart+off; 296 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&numPoints,&points);CHKERRQ(ierr); 297 } else { 298 ierr = PetscSectionGetOffset(csec,p,&off);CHKERRQ(ierr); 299 ierr = PetscSectionGetDof(csec,p,&dof);CHKERRQ(ierr); 300 for (q=0;q<dof/sdim;q++) vids[q] = off/sdim + q; 301 } 302 *nv = q; 303 PetscFunctionReturn(0); 304 } 305 306 static PetscErrorCode GLVisCreateFE(PetscFE femIn,char name[32],PetscFE *fem) 307 { 308 DM K; 309 PetscSpace P; 310 PetscDualSpace Q; 311 PetscQuadrature q,fq; 312 PetscInt dim,deg,dof; 313 DMPolytopeType ptype; 314 PetscBool isSimplex,isTensor; 315 PetscBool continuity = PETSC_FALSE; 316 PetscDTNodeType nodeType = PETSCDTNODES_GAUSSJACOBI; 317 PetscBool endpoint = PETSC_TRUE; 318 MPI_Comm comm; 319 PetscErrorCode ierr; 320 321 PetscFunctionBegin; 322 comm = PetscObjectComm((PetscObject)femIn); 323 ierr = PetscFEGetBasisSpace(femIn,&P);CHKERRQ(ierr); 324 ierr = PetscFEGetDualSpace(femIn,&Q);CHKERRQ(ierr); 325 ierr = PetscDualSpaceGetDM(Q,&K);CHKERRQ(ierr); 326 ierr = DMGetDimension(K,&dim);CHKERRQ(ierr); 327 ierr = PetscSpaceGetDegree(P,°,NULL);CHKERRQ(ierr); 328 ierr = PetscSpaceGetNumComponents(P,&dof);CHKERRQ(ierr); 329 ierr = DMPlexGetCellType(K,0,&ptype);CHKERRQ(ierr); 330 switch (ptype) { 331 case DM_POLYTOPE_QUADRILATERAL: 332 case DM_POLYTOPE_HEXAHEDRON: 333 isSimplex = PETSC_FALSE; break; 334 default: 335 isSimplex = PETSC_TRUE; break; 336 } 337 isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 338 /* Create space */ 339 ierr = PetscSpaceCreate(comm,&P);CHKERRQ(ierr); 340 ierr = PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 341 ierr = PetscSpacePolynomialSetTensor(P,isTensor);CHKERRQ(ierr); 342 ierr = PetscSpaceSetNumComponents(P,dof);CHKERRQ(ierr); 343 ierr = PetscSpaceSetNumVariables(P,dim);CHKERRQ(ierr); 344 ierr = PetscSpaceSetDegree(P,deg,PETSC_DETERMINE);CHKERRQ(ierr); 345 ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 346 /* Create dual space */ 347 ierr = PetscDualSpaceCreate(comm,&Q);CHKERRQ(ierr); 348 ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 349 ierr = PetscDualSpaceLagrangeSetTensor(Q,isTensor);CHKERRQ(ierr); 350 ierr = PetscDualSpaceLagrangeSetContinuity(Q,continuity);CHKERRQ(ierr); 351 ierr = PetscDualSpaceLagrangeSetNodeType(Q,nodeType,endpoint,0);CHKERRQ(ierr); 352 ierr = PetscDualSpaceSetNumComponents(Q,dof);CHKERRQ(ierr); 353 ierr = PetscDualSpaceSetOrder(Q,deg);CHKERRQ(ierr); 354 ierr = PetscDualSpaceCreateReferenceCell(Q,dim,isSimplex,&K);CHKERRQ(ierr); 355 ierr = PetscDualSpaceSetDM(Q,K);CHKERRQ(ierr); 356 ierr = DMDestroy(&K);CHKERRQ(ierr); 357 ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 358 /* Create quadrature */ 359 if (isSimplex) { 360 ierr = PetscDTStroudConicalQuadrature(dim, 1,deg+1,-1,+1,&q);CHKERRQ(ierr); 361 ierr = PetscDTStroudConicalQuadrature(dim-1,1,deg+1,-1,+1,&fq);CHKERRQ(ierr); 362 } else { 363 ierr = PetscDTGaussTensorQuadrature(dim, 1,deg+1,-1,+1,&q);CHKERRQ(ierr); 364 ierr = PetscDTGaussTensorQuadrature(dim-1,1,deg+1,-1,+1,&fq);CHKERRQ(ierr); 365 } 366 /* Create finite element */ 367 ierr = PetscFECreate(comm,fem);CHKERRQ(ierr); 368 ierr = PetscSNPrintf(name,32,"L2_T1_%DD_P%D",dim,deg);CHKERRQ(ierr); 369 ierr = PetscObjectSetName((PetscObject)*fem,name);CHKERRQ(ierr); 370 ierr = PetscFESetType(*fem,PETSCFEBASIC);CHKERRQ(ierr); 371 ierr = PetscFESetNumComponents(*fem,dof);CHKERRQ(ierr); 372 ierr = PetscFESetBasisSpace(*fem,P);CHKERRQ(ierr); 373 ierr = PetscFESetDualSpace(*fem,Q);CHKERRQ(ierr); 374 ierr = PetscFESetQuadrature(*fem,q);CHKERRQ(ierr); 375 ierr = PetscFESetFaceQuadrature(*fem,fq);CHKERRQ(ierr); 376 ierr = PetscFESetUp(*fem);CHKERRQ(ierr); 377 /* Cleanup */ 378 ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 379 ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 380 ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr); 381 ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr); 382 PetscFunctionReturn(0); 383 } 384 385 /* 386 ASCII visualization/dump: full support for simplices and tensor product cells. It supports AMR 387 Higher order meshes are also supported 388 */ 389 static PetscErrorCode DMPlexView_GLVis_ASCII(DM dm, PetscViewer viewer) 390 { 391 DMLabel label; 392 PetscSection coordSection,parentSection; 393 Vec coordinates,hovec; 394 const PetscScalar *array; 395 PetscInt bf,p,sdim,dim,depth,novl,minl; 396 PetscInt cStart,cEnd,vStart,vEnd,nvert; 397 PetscMPIInt size; 398 PetscBool localized,isascii; 399 PetscBool enable_mfem,enable_boundary,enable_ncmesh,view_ovl = PETSC_FALSE; 400 PetscBT pown,vown; 401 PetscErrorCode ierr; 402 PetscContainer glvis_container; 403 PetscBool cellvertex = PETSC_FALSE, periodic, enabled = PETSC_TRUE; 404 PetscBool enable_emark,enable_bmark; 405 const char *fmt; 406 char emark[64] = "",bmark[64] = ""; 407 408 PetscFunctionBegin; 409 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 410 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 411 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 412 if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII"); 413 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr); 414 if (size > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); 415 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 416 417 /* get container: determines if a process visualizes is portion of the data or not */ 418 ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr); 419 if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); 420 { 421 PetscViewerGLVisInfo glvis_info; 422 ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr); 423 enabled = glvis_info->enabled; 424 fmt = glvis_info->fmt; 425 } 426 427 /* Users can attach a coordinate vector to the DM in case they have a higher-order mesh 428 DMPlex does not currently support HO meshes, so there's no API for this */ 429 ierr = PetscObjectQuery((PetscObject)dm,"_glvis_mesh_coords",(PetscObject*)&hovec);CHKERRQ(ierr); 430 ierr = PetscObjectReference((PetscObject)hovec);CHKERRQ(ierr); 431 if (!hovec) { 432 DM cdm; 433 PetscFE disc; 434 PetscClassId classid; 435 436 ierr = DMGetCoordinateDM(dm,&cdm);CHKERRQ(ierr); 437 ierr = DMGetField(cdm,0,NULL,(PetscObject*)&disc);CHKERRQ(ierr); 438 ierr = PetscObjectGetClassId((PetscObject)disc,&classid);CHKERRQ(ierr); 439 if (classid == PETSCFE_CLASSID) { 440 DM hocdm; 441 PetscFE hodisc; 442 Vec vec; 443 Mat mat; 444 char name[32],fec_type[64]; 445 446 ierr = GLVisCreateFE(disc,name,&hodisc);CHKERRQ(ierr); 447 ierr = DMClone(cdm,&hocdm);CHKERRQ(ierr); 448 ierr = DMSetField(hocdm,0,NULL,(PetscObject)hodisc);CHKERRQ(ierr); 449 ierr = PetscFEDestroy(&hodisc);CHKERRQ(ierr); 450 ierr = DMCreateDS(hocdm);CHKERRQ(ierr); 451 452 ierr = DMGetCoordinates(dm,&vec);CHKERRQ(ierr); 453 ierr = DMCreateGlobalVector(hocdm,&hovec);CHKERRQ(ierr); 454 ierr = DMCreateInterpolation(cdm,hocdm,&mat,NULL);CHKERRQ(ierr); 455 ierr = MatInterpolate(mat,vec,hovec);CHKERRQ(ierr); 456 ierr = MatDestroy(&mat);CHKERRQ(ierr); 457 ierr = DMDestroy(&hocdm);CHKERRQ(ierr); 458 459 ierr = PetscSNPrintf(fec_type,sizeof(fec_type),"FiniteElementCollection: %s", name);CHKERRQ(ierr); 460 ierr = PetscObjectSetName((PetscObject)hovec,fec_type);CHKERRQ(ierr); 461 } 462 } 463 464 ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 465 ierr = DMPlexGetGhostCellStratum(dm,&p,NULL);CHKERRQ(ierr); 466 if (p >= 0) cEnd = p; 467 ierr = DMPlexGetDepthStratum(dm,0,&vStart,&vEnd);CHKERRQ(ierr); 468 ierr = DMGetPeriodicity(dm,&periodic,NULL,NULL,NULL);CHKERRQ(ierr); 469 ierr = DMGetCoordinatesLocalized(dm,&localized);CHKERRQ(ierr); 470 if (periodic && !localized && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Coordinates need to be localized"); 471 ierr = DMGetCoordinateSection(dm,&coordSection);CHKERRQ(ierr); 472 ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 473 ierr = DMGetCoordinatesLocal(dm,&coordinates);CHKERRQ(ierr); 474 if (!coordinates && !hovec) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Missing local coordinates vector"); 475 476 /* 477 a couple of sections of the mesh specification are disabled 478 - 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) 479 - vertex_parents: used for non-conforming meshes only when we want to use MFEM as a discretization package 480 and be able to derefine the mesh (MFEM does not currently have to ability to read ncmeshes in parallel) 481 */ 482 enable_boundary = PETSC_FALSE; 483 enable_ncmesh = PETSC_FALSE; 484 enable_mfem = PETSC_FALSE; 485 enable_emark = PETSC_FALSE; 486 enable_bmark = PETSC_FALSE; 487 /* I'm tired of problems with negative values in the markers, disable them */ 488 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)dm),((PetscObject)dm)->prefix,"GLVis PetscViewer DMPlex Options","PetscViewer");CHKERRQ(ierr); 489 ierr = PetscOptionsBool("-viewer_glvis_dm_plex_enable_boundary","Enable boundary section in mesh representation",NULL,enable_boundary,&enable_boundary,NULL);CHKERRQ(ierr); 490 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); 491 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); 492 ierr = PetscOptionsBool("-viewer_glvis_dm_plex_overlap","Include overlap region in local meshes",NULL,view_ovl,&view_ovl,NULL);CHKERRQ(ierr); 493 ierr = PetscOptionsString("-viewer_glvis_dm_plex_emarker","String for the material id label",NULL,emark,emark,sizeof(emark),&enable_emark);CHKERRQ(ierr); 494 ierr = PetscOptionsString("-viewer_glvis_dm_plex_bmarker","String for the boundary id label",NULL,bmark,bmark,sizeof(bmark),&enable_bmark);CHKERRQ(ierr); 495 ierr = PetscOptionsEnd();CHKERRQ(ierr); 496 if (enable_bmark) enable_boundary = PETSC_TRUE; 497 498 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 499 if (enable_ncmesh && size > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not supported in parallel"); 500 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 501 if (enable_boundary && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " 502 "Alternatively, run with -viewer_glvis_dm_plex_enable_boundary 0"); 503 if (enable_ncmesh && depth >= 0 && dim != depth) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mesh must be interpolated. " 504 "Alternatively, run with -viewer_glvis_dm_plex_enable_ncmesh 0"); 505 if (depth >=0 && dim != depth) { /* not interpolated, it assumes cell-vertex mesh */ 506 if (depth != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported depth %D. You should interpolate the mesh first",depth); 507 cellvertex = PETSC_TRUE; 508 } 509 510 /* Identify possible cells in the overlap */ 511 novl = 0; 512 pown = NULL; 513 if (size > 1) { 514 IS globalNum = NULL; 515 const PetscInt *gNum; 516 PetscBool ovl = PETSC_FALSE; 517 518 ierr = PetscObjectQuery((PetscObject)dm,"_glvis_plex_gnum",(PetscObject*)&globalNum);CHKERRQ(ierr); 519 if (!globalNum) { 520 if (view_ovl) { 521 ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),cEnd-cStart,0,1,&globalNum);CHKERRQ(ierr); 522 } else { 523 ierr = DMPlexCreateCellNumbering_Internal(dm,PETSC_TRUE,&globalNum);CHKERRQ(ierr); 524 } 525 ierr = PetscObjectCompose((PetscObject)dm,"_glvis_plex_gnum",(PetscObject)globalNum);CHKERRQ(ierr); 526 ierr = PetscObjectDereference((PetscObject)globalNum);CHKERRQ(ierr); 527 } 528 ierr = ISGetIndices(globalNum,&gNum);CHKERRQ(ierr); 529 for (p=cStart; p<cEnd; p++) { 530 if (gNum[p-cStart] < 0) { 531 ovl = PETSC_TRUE; 532 novl++; 533 } 534 } 535 if (ovl) { 536 /* it may happen that pown get not destroyed, if the user closes the window while this function is running. 537 TODO: garbage collector? attach pown to dm? */ 538 ierr = PetscBTCreate(cEnd-cStart,&pown);CHKERRQ(ierr); 539 for (p=cStart; p<cEnd; p++) { 540 if (gNum[p-cStart] < 0) continue; 541 else { 542 ierr = PetscBTSet(pown,p-cStart);CHKERRQ(ierr); 543 } 544 } 545 } 546 ierr = ISRestoreIndices(globalNum,&gNum);CHKERRQ(ierr); 547 } 548 549 /* return if this process is disabled */ 550 if (!enabled) { 551 ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 552 ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 553 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 554 ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 555 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 556 ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 557 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 558 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 559 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 560 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 561 ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 562 ierr = VecDestroy(&hovec);CHKERRQ(ierr); 563 PetscFunctionReturn(0); 564 } 565 566 if (enable_mfem) { 567 if (periodic && !hovec) { /* we need to generate a vector of L2 coordinates, as this is how MFEM handles periodic meshes */ 568 PetscInt vpc = 0; 569 char fec[64]; 570 PetscInt vids[8] = {0,1,2,3,4,5,6,7}; 571 PetscInt hexv[8] = {0,1,3,2,4,5,7,6}, tetv[4] = {0,1,2,3}; 572 PetscInt quadv[8] = {0,1,3,2}, triv[3] = {0,1,2}; 573 PetscInt *dof = NULL; 574 PetscScalar *array,*ptr; 575 576 ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: L2_T1_%DD_P1",dim);CHKERRQ(ierr); 577 if (cEnd-cStart) { 578 PetscInt fpc; 579 580 ierr = DMPlexGetConeSize(dm,cStart,&fpc);CHKERRQ(ierr); 581 switch(dim) { 582 case 1: 583 vpc = 2; 584 dof = hexv; 585 break; 586 case 2: 587 switch (fpc) { 588 case 3: 589 vpc = 3; 590 dof = triv; 591 break; 592 case 4: 593 vpc = 4; 594 dof = quadv; 595 break; 596 default: 597 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 598 } 599 break; 600 case 3: 601 switch (fpc) { 602 case 4: /* TODO: still need to understand L2 ordering for tets */ 603 vpc = 4; 604 dof = tetv; 605 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled tethraedral case"); 606 case 6: 607 if (cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: vertices per cell %D",fpc); 608 vpc = 8; 609 dof = hexv; 610 break; 611 case 8: 612 if (!cellvertex) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 613 vpc = 8; 614 dof = hexv; 615 break; 616 default: 617 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unhandled case: faces per cell %D",fpc); 618 } 619 break; 620 default: 621 SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unhandled dim"); 622 } 623 ierr = DMPlexReorderCell(dm,cStart,vids);CHKERRQ(ierr); 624 } 625 if (!dof) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing dofs"); 626 ierr = VecCreateSeq(PETSC_COMM_SELF,(cEnd-cStart-novl)*vpc*sdim,&hovec);CHKERRQ(ierr); 627 ierr = PetscObjectSetName((PetscObject)hovec,fec);CHKERRQ(ierr); 628 ierr = VecGetArray(hovec,&array);CHKERRQ(ierr); 629 ptr = array; 630 for (p=cStart;p<cEnd;p++) { 631 PetscInt csize,v,d; 632 PetscScalar *vals = NULL; 633 634 if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 635 ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 636 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); 637 for (v=0;v<vpc;v++) { 638 for (d=0;d<sdim;d++) { 639 ptr[sdim*dof[v]+d] = vals[sdim*vids[v]+d]; 640 } 641 } 642 ptr += vpc*sdim; 643 ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 644 } 645 ierr = VecRestoreArray(hovec,&array);CHKERRQ(ierr); 646 } 647 } 648 /* if we have high-order coordinates in 3D, we need to specify the boundary */ 649 if (hovec && dim == 3) enable_boundary = PETSC_TRUE; 650 651 /* header */ 652 ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 653 654 /* topological dimension */ 655 ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 656 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 657 658 /* elements */ 659 minl = 1; 660 label = NULL; 661 if (enable_emark) { 662 PetscInt lminl = PETSC_MAX_INT; 663 664 ierr = DMGetLabel(dm,emark,&label);CHKERRQ(ierr); 665 if (label) { 666 IS vals; 667 PetscInt ldef; 668 669 ierr = DMLabelGetDefaultValue(label,&ldef);CHKERRQ(ierr); 670 ierr = DMLabelGetValueIS(label,&vals);CHKERRQ(ierr); 671 ierr = ISGetMinMax(vals,&lminl,NULL);CHKERRQ(ierr); 672 ierr = ISDestroy(&vals);CHKERRQ(ierr); 673 lminl = PetscMin(ldef,lminl); 674 } 675 ierr = MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 676 if (minl == PETSC_MAX_INT) minl = 1; 677 } 678 ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 679 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",cEnd-cStart-novl);CHKERRQ(ierr); 680 for (p=cStart;p<cEnd;p++) { 681 PetscInt vids[8]; 682 PetscInt i,nv = 0,cid = -1,mid = 1; 683 684 if (PetscUnlikely(pown && !PetscBTLookup(pown,p-cStart))) continue; 685 ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);CHKERRQ(ierr); 686 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,p,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 687 ierr = DMPlexReorderCell(dm,p,vids);CHKERRQ(ierr); 688 ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 689 for (i=0;i<nv;i++) { 690 ierr = PetscViewerASCIIPrintf(viewer," %D",vids[i]);CHKERRQ(ierr); 691 } 692 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 693 } 694 695 /* boundary */ 696 ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 697 if (!enable_boundary) { 698 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 699 } else { 700 DMLabel perLabel; 701 PetscBT bfaces; 702 PetscInt fStart,fEnd,*fcells; 703 704 ierr = DMPlexGetHeightStratum(dm,1,&fStart,&fEnd);CHKERRQ(ierr); 705 ierr = PetscBTCreate(fEnd-fStart,&bfaces);CHKERRQ(ierr); 706 ierr = DMPlexGetMaxSizes(dm,NULL,&p);CHKERRQ(ierr); 707 ierr = PetscMalloc1(p,&fcells);CHKERRQ(ierr); 708 ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 709 if (!perLabel && localized) { /* this periodic cut can be moved up to DMPlex setup */ 710 ierr = DMCreateLabel(dm,"glvis_periodic_cut");CHKERRQ(ierr); 711 ierr = DMGetLabel(dm,"glvis_periodic_cut",&perLabel);CHKERRQ(ierr); 712 ierr = DMLabelSetDefaultValue(perLabel,1);CHKERRQ(ierr); 713 for (p=cStart;p<cEnd;p++) { 714 DMPolytopeType cellType; 715 PetscInt dof; 716 717 ierr = DMPlexGetCellType(dm,p,&cellType);CHKERRQ(ierr); 718 ierr = PetscSectionGetDof(coordSection,p,&dof);CHKERRQ(ierr); 719 if (dof) { 720 PetscInt uvpc, v,csize,cellClosureSize,*cellClosure = NULL,*vidxs = NULL; 721 PetscScalar *vals = NULL; 722 723 uvpc = DMPolytopeTypeGetNumVertices(cellType); 724 if (dof%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Incompatible number of cell dofs %D and space dimension %D",dof,sdim); 725 ierr = DMPlexVecGetClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 726 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 727 for (v=0;v<cellClosureSize;v++) 728 if (cellClosure[2*v] >= vStart && cellClosure[2*v] < vEnd) { 729 vidxs = cellClosure + 2*v; 730 break; 731 } 732 if (!vidxs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing vertices"); 733 for (v=0;v<uvpc;v++) { 734 PetscInt s; 735 736 for (s=0;s<sdim;s++) { 737 if (PetscAbsScalar(vals[v*sdim+s]-vals[v*sdim+s+uvpc*sdim])>PETSC_MACHINE_EPSILON) { 738 ierr = DMLabelSetValue(perLabel,vidxs[2*v],2);CHKERRQ(ierr); 739 } 740 } 741 } 742 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&cellClosureSize,&cellClosure);CHKERRQ(ierr); 743 ierr = DMPlexVecRestoreClosure(dm,coordSection,coordinates,p,&csize,&vals);CHKERRQ(ierr); 744 } 745 } 746 if (dim > 1) { 747 PetscInt eEnd,eStart; 748 749 ierr = DMPlexGetDepthStratum(dm,1,&eStart,&eEnd);CHKERRQ(ierr); 750 for (p=eStart;p<eEnd;p++) { 751 const PetscInt *cone; 752 PetscInt coneSize,i; 753 PetscBool ispe = PETSC_TRUE; 754 755 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 756 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 757 for (i=0;i<coneSize;i++) { 758 PetscInt v; 759 760 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 761 ispe = (PetscBool)(ispe && (v==2)); 762 } 763 if (ispe && coneSize) { 764 PetscInt ch, numChildren; 765 const PetscInt *children; 766 767 ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 768 ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 769 for (ch = 0; ch < numChildren; ch++) { 770 ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr); 771 } 772 } 773 } 774 if (dim > 2) { 775 for (p=fStart;p<fEnd;p++) { 776 const PetscInt *cone; 777 PetscInt coneSize,i; 778 PetscBool ispe = PETSC_TRUE; 779 780 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 781 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 782 for (i=0;i<coneSize;i++) { 783 PetscInt v; 784 785 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 786 ispe = (PetscBool)(ispe && (v==2)); 787 } 788 if (ispe && coneSize) { 789 PetscInt ch, numChildren; 790 const PetscInt *children; 791 792 ierr = DMLabelSetValue(perLabel,p,2);CHKERRQ(ierr); 793 ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 794 for (ch = 0; ch < numChildren; ch++) { 795 ierr = DMLabelSetValue(perLabel,children[ch],2);CHKERRQ(ierr); 796 } 797 } 798 } 799 } 800 } 801 } 802 for (p=fStart;p<fEnd;p++) { 803 const PetscInt *support; 804 PetscInt supportSize; 805 PetscBool isbf = PETSC_FALSE; 806 807 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 808 if (pown) { 809 PetscBool has_owned = PETSC_FALSE, has_ghost = PETSC_FALSE; 810 PetscInt i; 811 812 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 813 for (i=0;i<supportSize;i++) { 814 if (PetscLikely(PetscBTLookup(pown,support[i]-cStart))) has_owned = PETSC_TRUE; 815 else has_ghost = PETSC_TRUE; 816 } 817 isbf = (PetscBool)((supportSize == 1 && has_owned) || (supportSize > 1 && has_owned && has_ghost)); 818 } else { 819 isbf = (PetscBool)(supportSize == 1); 820 } 821 if (!isbf && perLabel) { 822 const PetscInt *cone; 823 PetscInt coneSize,i; 824 825 ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 826 ierr = DMPlexGetConeSize(dm,p,&coneSize);CHKERRQ(ierr); 827 isbf = PETSC_TRUE; 828 for (i=0;i<coneSize;i++) { 829 PetscInt v,d; 830 831 ierr = DMLabelGetValue(perLabel,cone[i],&v);CHKERRQ(ierr); 832 ierr = DMLabelGetDefaultValue(perLabel,&d);CHKERRQ(ierr); 833 isbf = (PetscBool)(isbf && v != d); 834 } 835 } 836 if (isbf) { 837 ierr = PetscBTSet(bfaces,p-fStart);CHKERRQ(ierr); 838 } 839 } 840 /* count boundary faces */ 841 for (p=fStart,bf=0;p<fEnd;p++) { 842 if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 843 const PetscInt *support; 844 PetscInt supportSize,c; 845 846 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 847 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 848 for (c=0;c<supportSize;c++) { 849 const PetscInt *cone; 850 PetscInt cell,cl,coneSize; 851 852 cell = support[c]; 853 if (pown && PetscUnlikely(!PetscBTLookup(pown,cell-cStart))) continue; 854 ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 855 ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 856 for (cl=0;cl<coneSize;cl++) { 857 if (cone[cl] == p) { 858 bf += 1; 859 break; 860 } 861 } 862 } 863 } 864 } 865 minl = 1; 866 label = NULL; 867 if (enable_bmark) { 868 PetscInt lminl = PETSC_MAX_INT; 869 870 ierr = DMGetLabel(dm,bmark,&label);CHKERRQ(ierr); 871 if (label) { 872 IS vals; 873 PetscInt ldef; 874 875 ierr = DMLabelGetDefaultValue(label,&ldef);CHKERRQ(ierr); 876 ierr = DMLabelGetValueIS(label,&vals);CHKERRQ(ierr); 877 ierr = ISGetMinMax(vals,&lminl,NULL);CHKERRQ(ierr); 878 ierr = ISDestroy(&vals);CHKERRQ(ierr); 879 lminl = PetscMin(ldef,lminl); 880 } 881 ierr = MPIU_Allreduce(&lminl,&minl,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 882 if (minl == PETSC_MAX_INT) minl = 1; 883 } 884 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",bf);CHKERRQ(ierr); 885 for (p=fStart;p<fEnd;p++) { 886 if (PetscUnlikely(PetscBTLookup(bfaces,p-fStart))) { 887 const PetscInt *support; 888 PetscInt supportSize,c,nc = 0; 889 890 ierr = DMPlexGetSupportSize(dm,p,&supportSize);CHKERRQ(ierr); 891 ierr = DMPlexGetSupport(dm,p,&support);CHKERRQ(ierr); 892 if (pown) { 893 for (c=0;c<supportSize;c++) { 894 if (PetscLikely(PetscBTLookup(pown,support[c]-cStart))) { 895 fcells[nc++] = support[c]; 896 } 897 } 898 } else for (c=0;c<supportSize;c++) fcells[nc++] = support[c]; 899 for (c=0;c<nc;c++) { 900 const DMPolytopeType *faceTypes; 901 DMPolytopeType cellType; 902 const PetscInt *faceSizes,*cone; 903 PetscInt vids[8],*faces,st,i,coneSize,cell,cl,nv,cid = -1,mid = -1; 904 905 cell = fcells[c]; 906 ierr = DMPlexGetCone(dm,cell,&cone);CHKERRQ(ierr); 907 ierr = DMPlexGetConeSize(dm,cell,&coneSize);CHKERRQ(ierr); 908 for (cl=0;cl<coneSize;cl++) 909 if (cone[cl] == p) 910 break; 911 if (cl == coneSize) continue; 912 913 /* face material id and type */ 914 ierr = DMPlexGetPointMFEMCellID_Internal(dm,label,minl,p,&mid,&cid);CHKERRQ(ierr); 915 ierr = PetscViewerASCIIPrintf(viewer,"%D %D",mid,cid);CHKERRQ(ierr); 916 /* vertex ids */ 917 ierr = DMPlexGetCellType(dm,cell,&cellType);CHKERRQ(ierr); 918 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,cell,(localized && !hovec) ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 919 ierr = DMPlexGetRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces);CHKERRQ(ierr); 920 st = 0; 921 for (i=0;i<cl;i++) st += faceSizes[i]; 922 ierr = DMPlexInvertCell(faceTypes[cl],faces + st);CHKERRQ(ierr); 923 for (i=0;i<faceSizes[cl];i++) { 924 ierr = PetscViewerASCIIPrintf(viewer," %d",faces[st+i]);CHKERRQ(ierr); 925 } 926 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 927 ierr = DMPlexRestoreRawFaces_Internal(dm,cellType,vids,NULL,&faceTypes,&faceSizes,(const PetscInt**)&faces);CHKERRQ(ierr); 928 bf -= 1; 929 } 930 } 931 } 932 if (bf) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Remaining boundary faces %D",bf); 933 ierr = PetscBTDestroy(&bfaces);CHKERRQ(ierr); 934 ierr = PetscFree(fcells);CHKERRQ(ierr); 935 } 936 937 /* mark owned vertices */ 938 vown = NULL; 939 if (pown) { 940 ierr = PetscBTCreate(vEnd-vStart,&vown);CHKERRQ(ierr); 941 for (p=cStart;p<cEnd;p++) { 942 PetscInt i,closureSize,*closure = NULL; 943 944 if (PetscUnlikely(!PetscBTLookup(pown,p-cStart))) continue; 945 ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 946 for (i=0;i<closureSize;i++) { 947 const PetscInt pp = closure[2*i]; 948 949 if (pp >= vStart && pp < vEnd) { 950 ierr = PetscBTSet(vown,pp-vStart);CHKERRQ(ierr); 951 } 952 } 953 ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 954 } 955 } 956 957 /* vertex_parents (Non-conforming meshes) */ 958 parentSection = NULL; 959 if (enable_ncmesh) { 960 ierr = DMPlexGetTree(dm,&parentSection,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 961 } 962 if (parentSection) { 963 PetscInt vp,gvp; 964 965 for (vp=0,p=vStart;p<vEnd;p++) { 966 DMLabel dlabel; 967 PetscInt parent,depth; 968 969 if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 970 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 971 ierr = DMLabelGetValue(dlabel,p,&depth);CHKERRQ(ierr); 972 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 973 if (parent != p) vp++; 974 } 975 ierr = MPIU_Allreduce(&vp,&gvp,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 976 if (gvp) { 977 PetscInt maxsupp; 978 PetscBool *skip = NULL; 979 980 ierr = PetscViewerASCIIPrintf(viewer,"\nvertex_parents\n");CHKERRQ(ierr); 981 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vp);CHKERRQ(ierr); 982 ierr = DMPlexGetMaxSizes(dm,NULL,&maxsupp);CHKERRQ(ierr); 983 ierr = PetscMalloc1(maxsupp,&skip);CHKERRQ(ierr); 984 for (p=vStart;p<vEnd;p++) { 985 DMLabel dlabel; 986 PetscInt parent; 987 988 if (PetscUnlikely(vown && !PetscBTLookup(vown,p-vStart))) continue; 989 ierr = DMPlexGetDepthLabel(dm,&dlabel);CHKERRQ(ierr); 990 ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 991 if (parent != p) { 992 PetscInt vids[8] = { -1, -1, -1, -1, -1, -1, -1, -1 }; /* silent overzealous clang static analyzer */ 993 PetscInt i,nv,ssize,n,numChildren,depth = -1; 994 const PetscInt *children; 995 996 ierr = DMPlexGetConeSize(dm,parent,&ssize);CHKERRQ(ierr); 997 switch (ssize) { 998 case 2: /* edge */ 999 nv = 0; 1000 ierr = DMPlexGetPointMFEMVertexIDs_Internal(dm,parent,localized ? coordSection : NULL,&nv,vids);CHKERRQ(ierr); 1001 ierr = PetscViewerASCIIPrintf(viewer,"%D",p-vStart);CHKERRQ(ierr); 1002 for (i=0;i<nv;i++) { 1003 ierr = PetscViewerASCIIPrintf(viewer," %D",vids[i]);CHKERRQ(ierr); 1004 } 1005 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1006 vp--; 1007 break; 1008 case 4: /* face */ 1009 ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr); 1010 for (n=0;n<numChildren;n++) { 1011 ierr = DMLabelGetValue(dlabel,children[n],&depth);CHKERRQ(ierr); 1012 if (!depth) { 1013 const PetscInt *hvsupp,*hesupp,*cone; 1014 PetscInt hvsuppSize,hesuppSize,coneSize; 1015 PetscInt hv = children[n],he = -1,f; 1016 1017 ierr = PetscArrayzero(skip,maxsupp);CHKERRQ(ierr); 1018 ierr = DMPlexGetSupportSize(dm,hv,&hvsuppSize);CHKERRQ(ierr); 1019 ierr = DMPlexGetSupport(dm,hv,&hvsupp);CHKERRQ(ierr); 1020 for (i=0;i<hvsuppSize;i++) { 1021 PetscInt ep; 1022 ierr = DMPlexGetTreeParent(dm,hvsupp[i],&ep,NULL);CHKERRQ(ierr); 1023 if (ep != hvsupp[i]) { 1024 he = hvsupp[i]; 1025 } else { 1026 skip[i] = PETSC_TRUE; 1027 } 1028 } 1029 if (he == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vertex %D support size %D: hanging edge not found",hv,hvsuppSize); 1030 ierr = DMPlexGetCone(dm,he,&cone);CHKERRQ(ierr); 1031 vids[0] = (cone[0] == hv) ? cone[1] : cone[0]; 1032 ierr = DMPlexGetSupportSize(dm,he,&hesuppSize);CHKERRQ(ierr); 1033 ierr = DMPlexGetSupport(dm,he,&hesupp);CHKERRQ(ierr); 1034 for (f=0;f<hesuppSize;f++) { 1035 PetscInt j; 1036 1037 ierr = DMPlexGetCone(dm,hesupp[f],&cone);CHKERRQ(ierr); 1038 ierr = DMPlexGetConeSize(dm,hesupp[f],&coneSize);CHKERRQ(ierr); 1039 for (j=0;j<coneSize;j++) { 1040 PetscInt k; 1041 for (k=0;k<hvsuppSize;k++) { 1042 if (hvsupp[k] == cone[j]) { 1043 skip[k] = PETSC_TRUE; 1044 break; 1045 } 1046 } 1047 } 1048 } 1049 for (i=0;i<hvsuppSize;i++) { 1050 if (!skip[i]) { 1051 ierr = DMPlexGetCone(dm,hvsupp[i],&cone);CHKERRQ(ierr); 1052 vids[1] = (cone[0] == hv) ? cone[1] : cone[0]; 1053 } 1054 } 1055 ierr = PetscViewerASCIIPrintf(viewer,"%D",hv-vStart);CHKERRQ(ierr); 1056 for (i=0;i<2;i++) { 1057 ierr = PetscViewerASCIIPrintf(viewer," %D",vids[i]-vStart);CHKERRQ(ierr); 1058 } 1059 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1060 vp--; 1061 } 1062 } 1063 break; 1064 default: 1065 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Don't know how to deal with support size %D",ssize); 1066 } 1067 } 1068 } 1069 ierr = PetscFree(skip);CHKERRQ(ierr); 1070 } 1071 if (vp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unexpected %D hanging vertices",vp); 1072 } 1073 ierr = PetscBTDestroy(&pown);CHKERRQ(ierr); 1074 ierr = PetscBTDestroy(&vown);CHKERRQ(ierr); 1075 1076 /* vertices */ 1077 if (hovec) { /* higher-order meshes */ 1078 const char *fec; 1079 PetscInt i,n,s; 1080 1081 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 1082 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",vEnd-vStart);CHKERRQ(ierr); 1083 ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr); 1084 ierr = PetscObjectGetName((PetscObject)hovec,&fec);CHKERRQ(ierr); 1085 ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr); 1086 ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fec);CHKERRQ(ierr); 1087 ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr); 1088 ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/ 1089 ierr = VecGetArrayRead(hovec,&array);CHKERRQ(ierr); 1090 ierr = VecGetLocalSize(hovec,&n);CHKERRQ(ierr); 1091 if (n%sdim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Size of local coordinate vector %D incompatible with space dimension %D",n,sdim); 1092 for (i=0;i<n/sdim;i++) { 1093 for (s=0;s<sdim;s++) { 1094 ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[i*sdim+s]));CHKERRQ(ierr); 1095 } 1096 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1097 } 1098 ierr = VecRestoreArrayRead(hovec,&array);CHKERRQ(ierr); 1099 } else { 1100 ierr = VecGetLocalSize(coordinates,&nvert);CHKERRQ(ierr); 1101 ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 1102 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nvert/sdim);CHKERRQ(ierr); 1103 ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 1104 ierr = VecGetArrayRead(coordinates,&array);CHKERRQ(ierr); 1105 for (p=0;p<nvert/sdim;p++) { 1106 PetscInt s; 1107 for (s=0;s<sdim;s++) { 1108 PetscReal v = PetscRealPart(array[p*sdim+s]); 1109 1110 ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscIsInfOrNanReal(v) ? 0.0 : v);CHKERRQ(ierr); 1111 } 1112 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1113 } 1114 ierr = VecRestoreArrayRead(coordinates,&array);CHKERRQ(ierr); 1115 } 1116 ierr = VecDestroy(&hovec);CHKERRQ(ierr); 1117 PetscFunctionReturn(0); 1118 } 1119 1120 PetscErrorCode DMPlexView_GLVis(DM dm, PetscViewer viewer) 1121 { 1122 PetscErrorCode ierr; 1123 PetscFunctionBegin; 1124 ierr = DMView_GLVis(dm,viewer,DMPlexView_GLVis_ASCII);CHKERRQ(ierr); 1125 PetscFunctionReturn(0); 1126 } 1127