/* * *************************************************************************** * Routine: Gem_countFaces * * Purpose: Count all of the faces without actually creating them. * * Notes: Keep track of the global face numbers in each element. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_countFaces(Gem *thee) { int i, j, faceCnt, face; # if defined(VG_ELEMENT) int face2; # endif SS *sm, *sm2; Vnm_tstart(70, "face count"); Vnm_print(0,"Gem_countFaces: counting faces.."); /* count the faces (including boundary faces) */ if (Gem_dim(thee) == 2) { faceCnt = 0; } else { faceCnt = 0; for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); /* now count all i-j (sm-sm2) faces where i<j */ for (face=0; face<Gem_dimVV(thee); face++) { sm2 = SS_nabor(sm,face); if (sm2 == VNULL) { # if defined(VG_ELEMENT) SS_setFaceNumber(sm,face,faceCnt); # endif faceCnt++; } else { j = SS_id(sm2); if (i < j) { # if defined(VG_ELEMENT) face2 = SS_sharedFaceLocNum(sm2,sm); SS_setFaceNumber(sm,face,faceCnt); SS_setFaceNumber(sm2,face2,faceCnt); # endif faceCnt++; } } } } } Vnm_print(0,"..done. [faces=<%d>]\n",faceCnt); Vnm_tstop(70, "face count"); /* set the face counter and return */ Gem_setNumVirtFF(thee, faceCnt); }
VPUBLIC int Vcsm_update(Vcsm *thee, SS **simps, int num) { /* Counters */ int isimp, jsimp, iatom, jatom, atomID, simpID; int nsimps, gotMem; /* Object info */ Vatom *atom; SS *simplex; double *position; /* Lists */ int *qParent, nqParent; int **sqmNew, *nsqmNew; int *affAtoms, nAffAtoms; int *dnqsm, *nqsmNew, **qsmNew; VASSERT(thee != VNULL); VASSERT(thee->initFlag); /* If we don't have enough memory to accommodate the new entries, * add more by doubling the existing amount */ isimp = thee->nsimp + num - 1; gotMem = 0; while (!gotMem) { if (isimp > thee->msimp) { isimp = 2 * isimp; thee->nsqm = Vmem_realloc(thee->vmem, thee->msimp, sizeof(int), (void **)&(thee->nsqm), isimp); VASSERT(thee->nsqm != VNULL); thee->sqm = Vmem_realloc(thee->vmem, thee->msimp, sizeof(int *), (void **)&(thee->sqm), isimp); VASSERT(thee->sqm != VNULL); thee->msimp = isimp; } else gotMem = 1; } /* Initialize the nsqm entires we just allocated */ for (isimp = thee->nsimp; isimp<thee->nsimp+num-1 ; isimp++) { thee->nsqm[isimp] = 0; } thee->nsimp = thee->nsimp + num - 1; /* There's a simple case to deal with: if simps[0] didn't have a * charge in the first place */ isimp = SS_id(simps[0]); if (thee->nsqm[isimp] == 0) { for (isimp=1; isimp<num; isimp++) { thee->nsqm[SS_id(simps[isimp])] = 0; } return 1; } /* The more complicated case has occured; the parent simplex had one or * more charges. First, generate the list of affected charges. */ isimp = SS_id(simps[0]); nqParent = thee->nsqm[isimp]; qParent = thee->sqm[isimp]; sqmNew = Vmem_malloc(thee->vmem, num, sizeof(int *)); VASSERT(sqmNew != VNULL); nsqmNew = Vmem_malloc(thee->vmem, num, sizeof(int)); VASSERT(nsqmNew != VNULL); for (isimp=0; isimp<num; isimp++) nsqmNew[isimp] = 0; /* Loop throught the affected atoms to determine how many atoms each * simplex will get. */ for (iatom=0; iatom<nqParent; iatom++) { atomID = qParent[iatom]; atom = Valist_getAtom(thee->alist, atomID); position = Vatom_getPosition(atom); nsimps = 0; jsimp = 0; for (isimp=0; isimp<num; isimp++) { simplex = simps[isimp]; if (Gem_pointInSimplex(thee->gm, simplex, position)) { nsqmNew[isimp]++; jsimp = 1; } } VASSERT(jsimp != 0); } /* Sanity check that we didn't lose any atoms... */ iatom = 0; for (isimp=0; isimp<num; isimp++) iatom += nsqmNew[isimp]; if (iatom < nqParent) { Vnm_print(2,"Vcsm_update: Lost %d (of %d) atoms!\n", nqParent - iatom, nqParent); VASSERT(0); } /* Allocate the storage */ for (isimp=0; isimp<num; isimp++) { if (nsqmNew[isimp] > 0) { sqmNew[isimp] = Vmem_malloc(thee->vmem, nsqmNew[isimp], sizeof(int)); VASSERT(sqmNew[isimp] != VNULL); } } /* Assign charges to simplices */ for (isimp=0; isimp<num; isimp++) { jsimp = 0; simplex = simps[isimp]; /* Loop over the atoms associated with the parent simplex */ for (iatom=0; iatom<nqParent; iatom++) { atomID = qParent[iatom]; atom = Valist_getAtom(thee->alist, atomID); position = Vatom_getPosition(atom); if (Gem_pointInSimplex(thee->gm, simplex, position)) { sqmNew[isimp][jsimp] = atomID; jsimp++; } } } /* Update the QSM map using the old and new SQM lists */ /* The affected atoms are those contained in the parent simplex; i.e. * thee->sqm[SS_id(simps[0])] */ affAtoms = thee->sqm[SS_id(simps[0])]; nAffAtoms = thee->nsqm[SS_id(simps[0])]; /* Each of these atoms will go somewhere else; i.e., the entries in * thee->qsm are never destroyed and thee->nqsm never decreases. * However, it is possible that a subdivision could cause an atom to be * shared by two child simplices. Here we record the change, if any, * in the number of simplices associated with each atom. */ dnqsm = Vmem_malloc(thee->vmem, nAffAtoms, sizeof(int)); VASSERT(dnqsm != VNULL); nqsmNew = Vmem_malloc(thee->vmem, nAffAtoms, sizeof(int)); VASSERT(nqsmNew != VNULL); qsmNew = Vmem_malloc(thee->vmem, nAffAtoms, sizeof(int*)); VASSERT(qsmNew != VNULL); for (iatom=0; iatom<nAffAtoms; iatom++) { dnqsm[iatom] = -1; atomID = affAtoms[iatom]; for (isimp=0; isimp<num; isimp++) { for (jatom=0; jatom<nsqmNew[isimp]; jatom++) { if (sqmNew[isimp][jatom] == atomID) dnqsm[iatom]++; } } VASSERT(dnqsm[iatom] > -1); } /* Setup the new entries in the array */ for (iatom=0;iatom<nAffAtoms; iatom++) { atomID = affAtoms[iatom]; qsmNew[iatom] = Vmem_malloc(thee->vmem, (dnqsm[iatom] + thee->nqsm[atomID]), sizeof(int)); nqsmNew[iatom] = 0; VASSERT(qsmNew[iatom] != VNULL); } /* Fill the new entries in the array */ /* First, do the modified entries */ for (isimp=0; isimp<num; isimp++) { simpID = SS_id(simps[isimp]); for (iatom=0; iatom<nsqmNew[isimp]; iatom++) { atomID = sqmNew[isimp][iatom]; for (jatom=0; jatom<nAffAtoms; jatom++) { if (atomID == affAtoms[jatom]) break; } if (jatom < nAffAtoms) { qsmNew[jatom][nqsmNew[jatom]] = simpID; nqsmNew[jatom]++; } } } /* Now do the unmodified entries */ for (iatom=0; iatom<nAffAtoms; iatom++) { atomID = affAtoms[iatom]; for (isimp=0; isimp<thee->nqsm[atomID]; isimp++) { for (jsimp=0; jsimp<num; jsimp++) { simpID = SS_id(simps[jsimp]); if (thee->qsm[atomID][isimp] == simpID) break; } if (jsimp == num) { qsmNew[iatom][nqsmNew[iatom]] = thee->qsm[atomID][isimp]; nqsmNew[iatom]++; } } } /* Replace the existing entries in the table. Do the QSM entires * first, since they require affAtoms = thee->sqm[simps[0]] */ for (iatom=0; iatom<nAffAtoms; iatom++) { atomID = affAtoms[iatom]; Vmem_free(thee->vmem, thee->nqsm[atomID], sizeof(int), (void **)&(thee->qsm[atomID])); thee->qsm[atomID] = qsmNew[iatom]; thee->nqsm[atomID] = nqsmNew[iatom]; } for (isimp=0; isimp<num; isimp++) { simpID = SS_id(simps[isimp]); if (thee->nsqm[simpID] > 0) Vmem_free(thee->vmem, thee->nsqm[simpID], sizeof(int), (void **)&(thee->sqm[simpID])); thee->sqm[simpID] = sqmNew[isimp]; thee->nsqm[simpID] = nsqmNew[isimp]; } Vmem_free(thee->vmem, num, sizeof(int *), (void **)&sqmNew); Vmem_free(thee->vmem, num, sizeof(int), (void **)&nsqmNew); Vmem_free(thee->vmem, nAffAtoms, sizeof(int *), (void **)&qsmNew); Vmem_free(thee->vmem, nAffAtoms, sizeof(int), (void **)&nqsmNew); Vmem_free(thee->vmem, nAffAtoms, sizeof(int), (void **)&dnqsm); return 1; }
/* * *************************************************************************** * Routine: Gem_contentChk * * Purpose: Print out contents of all geometry structures. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_contentChk(Gem *thee) { int i, j, k, ii; SS *sm; EE *eg; VV *vx; Vnm_print(0,"Gem_contentChk: printing data.\n"); for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); Vnm_print(0,"SIMPLEX [%d]\n", SS_id(sm)); Vnm_print(0," g.bits =<%d>\n",sm->g.bits); Vnm_print(0," g.uid =<%d>\n",sm->g.uid); Vnm_print(0," d.tags =<%d>\n",sm->d.tags); Vnm_print(0," d.faces=<%d>\n",sm->d.faces); Vnm_print(0," d.vPtr ="); for (ii=0; ii<4; ii++) { j = -1; if (sm->d.vPtr[ii] != VNULL) j = VV_id(sm->d.vPtr[ii]); Vnm_print(0," <%d>",j); } Vnm_print(0,"\n"); Vnm_print(0," d.sPtr ="); for (ii=0; ii<4; ii++) { j = -1; if (sm->d.sPtr[ii] != VNULL) j = SS_id(sm->d.sPtr[ii]); Vnm_print(0," <%d>",j); } Vnm_print(0,"\n"); # if defined(VG_ELEMENT) Vnm_print(0," d.fNum ="); for (ii=0; ii<4; ii++) { j = sm->d.fNum[ii]; Vnm_print(0," <%d>",j); } Vnm_print(0,"\n"); Vnm_print(0," d.eNum ="); for (ii=0; ii<6; ii++) { j = sm->d.eNum[ii]; Vnm_print(0," <%d>",j); } Vnm_print(0,"\n"); # endif } for (i=0; i<Gem_numEE(thee); i++) { eg = Gem_EE(thee,i); Vnm_print(0,"EDGE [%d]\n", EE_id(eg)); Vnm_print(0," g.bits=<%d>\n",eg->g.bits); Vnm_print(0," g.uid =<%d>\n",eg->g.uid); k = 0; if (eg->d.midPtr != VNULL) k = VV_id(eg->d.midPtr); Vnm_print(0," d.midPtr=<%d>\n", k); Vnm_print(0," d.vPtr ="); for (ii=0; ii<2; ii++) { j = 0; if (eg->d.vPtr[ii] != 0) j = VV_id(eg->d.vPtr[ii]); Vnm_print(0," <%d>",j); } Vnm_print(0,"\n"); Vnm_print(0," d.ePtr ="); for (ii=0; ii<2; ii++) { j = 0; if (eg->d.ePtr[ii] != 0) j = EE_id(eg->d.ePtr[ii]); Vnm_print(0," <%d>",j); } Vnm_print(0,"\n"); } for (i=0; i<Gem_numVV(thee); i++) { vx = Gem_VV(thee,i); Vnm_print(0,"VERTEX [%d]\n", VV_id(vx)); Vnm_print(0," g.bits=<%d>\n",vx->g.bits); Vnm_print(0," g.uid =<%d>\n",vx->g.uid); j = 0; k = 0; if (vx->d.ePtr != VNULL) j = EE_id(vx->d.ePtr); if (vx->d.sPtr != VNULL) k = SS_id(vx->d.sPtr); Vnm_print(0," d.ePtr=<%d>\n", j); Vnm_print(0," d.sPtr=<%d>\n", k); Vnm_print(0," d.x ="); for (ii=0; ii<3; ii++) Vnm_print(0," <%e>", vx->d.x[ii]); Vnm_print(0,"\n"); } }
/* * *************************************************************************** * Routine: Aprx_partSpect * * Purpose: Partition the domain using spectral bisection. * * Notes: We solve the following generalized eigenvalue problem: * * Ax = lambda Bx * * for second smallest eigenpair. We then return the eigenvector * from the pair. We make this happen by turning it into a * regular eigenvalue problem: * * B^{-1/2} A B^{-1/2} ( B^{1/2} x ) = lambda ( B^{1/2} x ) * * or rather * * C y = lambda y, where C=B^{-1/2}AB^{-1/2}, y=B^{1/2}x. * * The matrix "B" is simply a diagonal matrix with a (positive) * error estimate for the element on the diagonal. Therefore, the * matrix B^{-1/2} is a well-defined positive diagonal matrix. * * We explicitly form the matrix C and then send it to the inverse * Rayleigh-quotient iteration to recover the second smallest * eigenpair. On return, we scale the eigenvector y by B^{-1/2} to * recover the actual eigenvector x = B^{-1/2} y. * * To handle the possibility that an element has zero error, in * which case the B matrix had a zero on the diagonal, we set the * corresponding entry in B^{-1/2} to be 1. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC int Aprx_partSpect(Aprx *thee, int pcolor, int numC, double *evec, simHelper *simH, int *ford, int *rord, int general) { int i, j, k, dim, itmax, litmax, key, flag; int numF, *IA, *JA; int numB, numR[MAXV]; MATsym psym[MAXV][MAXV]; MATmirror pmirror[MAXV][MAXV]; MATformat pfrmt[MAXV][MAXV]; int numO[MAXV][MAXV], *IJA[MAXV][MAXV]; double lambda, normal, etol, letol, value; SS *sm, *sm0; VV *v[4]; Bmat *A; Bvec *u, *B2, *B2inv; Vnm_print(0,"Aprx_partSpect: [pc=%d] partitioning:\n", pcolor); /* dimensions */ dim = Gem_dim(thee->gm); /* go through the elements and enumerate elements and faces */ numF=0; for (i=0; i<numC; i++) { sm = Gem_SS(thee->gm,ford[i]); for (j=0;j<(int)SS_dimVV(sm);j++) { /* get the vertex pointers for this face */ for (k=0; k<dim; k++) { v[k] = SS_vertex( sm, vmapF[j][k] ); } /* do we already know our nabor */ if (dim == 2) { /* find the unique face nabor (if not on boundary) */ sm0 = VV_commonSimplex2(v[0],v[1],sm); } else { /* (dim == 3) */ /* find the unique face nabor (if not on boundary) */ sm0 = VV_commonSimplex3(v[0],v[1],v[2],sm); } /* okay we found a nabor */ if (sm0 != VNULL) { if ((int)SS_chart(sm0) == pcolor) { simH[i].diag++; k = rord[ SS_id(sm0) ]; if (k > i) { simH[i].faceId[j] = k; numF++; } } } } } /* el; loop over elements */ /* sort the rows and build matrix integer structures */ IJA[0][0] = Vmem_malloc( thee->vmem, numC+1+numF, sizeof(int) ); IA = IJA[0][0]; JA = IA + numC + 1; k = 0; IA[0] = 0; for (i=0; i<numC; i++) { Vnm_qsort(simH[i].faceId, 4); for (j=0; j<4; j++) { if (simH[i].faceId[j] > i) { JA[k] = simH[i].faceId[j]; k++; } } IA[i+1] = k; } VASSERT( k == numF ); /* create the real matrix object (creating space for matrix entries) */ numB = 1; numR[0] = numC; numO[0][0] = numF; psym[0][0] = IS_SYM; /* symmetric */ pmirror[0][0] = ISNOT_MIRROR; /* really exists */ pfrmt[0][0] = DRC_FORMAT; /* YSMP-bank */ A = Bmat_ctor( thee->vmem, "A", numB, numR, numR, pmirror ); Bmat_initStructure( A, pfrmt, psym, numO, IJA ); /* create the eigenvector */ u = Bvec_ctor( thee->vmem, "u", numB, numR ); /* create the scaling matrix */ B2 = Bvec_ctor( thee->vmem, "B2", numB, numR ); B2inv = Bvec_ctor( thee->vmem, "B2inv", numB, numR ); for (i=0; i<numC; i++) { if ( general ) { if ( simH[i].error > 0. ) { Bvec_setB( B2, 0, i, VSQRT(simH[i].error) ); Bvec_setB( B2inv, 0, i, 1./VSQRT(simH[i].error) ); } else if ( simH[i].error == 0. ) { Bvec_setB( B2, 0, i, 1. ); Bvec_setB( B2inv, 0, i, 1. ); } else { VASSERT(0); } } else { Bvec_setB( B2, 0, i, 1. ); Bvec_setB( B2inv, 0, i, 1. ); } } /* now build the scaled adjacency matrix entries */ k = 0; for (i=0; i<numC; i++) { value = (double)simH[i].diag * Bvec_valB( B2inv, 0, i ) * Bvec_valB( B2inv, 0, i ); Bmat_set( A, 0, 0, i, i, value ); for (j=0; j<4; j++) { if (simH[i].faceId[j] > i) { value = -1. * Bvec_valB( B2inv, 0, i ) * Bvec_valB( B2inv, 0, simH[i].faceId[j] ); Bmat_set( A, 0, 0, i, simH[i].faceId[j], value ); k++; } } } VASSERT( k == numO[0][0] ); /* the initial approximation */ for (i=0; i<numC; i++) { Bvec_setB( u, 0, i, evec[i] ); } /* print out matrix for testing */ /* Bmat_printSp(A, "lap.m"); */ /* finally, get eigenvector #2 (this is the costly part...) */ litmax = 500; letol = 1.0e-3; lambda = 0.; key = 0; flag = 0; itmax = 50; etol = 1.0e-4; Bvec_eig(u, A, litmax, letol, &lambda, key, flag, itmax, etol); /* re-scale the final approximation */ normal = 0; for (i=0; i<numC; i++) { evec[i] = Bvec_valB( B2inv, 0, i ) * Bvec_valB( u, 0, i ); normal += (evec[i]*evec[i]); } normal = VSQRT( normal ); /* normalize the final result */ for (i=0; i<numC; i++) { evec[i] = evec[i] / normal; } /* destroy the adjacency matrix and eigenvector */ Bmat_dtor( &A ); /* this frees our earlier IJA malloc! */ Bvec_dtor( &B2 ); Bvec_dtor( &B2inv ); Bvec_dtor( &u ); return 0; }
/* * *************************************************************************** * Routine: Gem_formChk * * Purpose: Check the self-consistency of the geometry datastructures. * * Notes: key==0 --> check: min (just vertices and simplices) * key==1 --> check: min + simplex ring * key==2 --> check: min + simplex ring + edge ring * key==3 --> check: min + simplex ring + edge ring + conform * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_formChk(Gem *thee, int key) { int i, j, k; int l, m, bndVert, bndFace; int edgeOrderProb, conformCount, fType[4]; double svol; VV *vx, *v[4]; SS *sm, *sm0, *sm1, *sm2; EE *eg; /* input check and some i/o */ VASSERT( (key >= 0) && (key <= 3) ); /* go through all vertices and check for consistency */ Vnm_tstart(80, "form check"); Vnm_print(0,"Gem_formChk: Topology check on: " "<%d> SS, <%d> EE, <%d> VV:\n", Gem_numSS(thee), Gem_numEE(thee), Gem_numVV(thee)); Vnm_print(0,"Gem_formChk: ..VV..\n"); bndVert = 0; for (i=0; i<Gem_numVV(thee); i++) { vx = Gem_VV(thee,i); if ( (i>0) && (i % VPRTKEY) == 0 ) Vnm_print(0,"[CV:%d]",i); VJMPERR1( !Vnm_sigInt() ); /* check basic data */ if ((int)VV_id(vx)!=i) Vnm_print(2,"Gem_formChk: VV_id BAD vtx <%d>\n", i); if (key == 0) { if (VV_firstSS(vx) != VNULL) Vnm_print(2,"Gem_formChk: firstS BAD vtx <%d>\n",i); if (VV_firstEE(vx) != VNULL) Vnm_print(2,"Gem_formChk: firstE BAD vtx <%d>\n",i); } else if (key == 1) { if (VV_firstSS(vx) == VNULL) Vnm_print(2,"Gem_formChk: firstS BAD vtx <%d>\n",i); if (VV_firstEE(vx) != VNULL) Vnm_print(2,"Gem_formChk: firstE BAD vtx <%d>\n",i); } else if ((key == 2) || (key == 3)) { if (VV_firstSS(vx) == VNULL) Vnm_print(2,"Gem_formChk: firstS BAD vtx <%d>\n",i); if ((VV_firstEE(vx) == VNULL) && (Gem_numEE(thee) > 0)) Vnm_print(2,"Gem_formChk: firstE BAD vtx <%d>\n",i); } else { VASSERT(0); } /* boundary data */ if ( VBOUNDARY( VV_type(vx) ) ) bndVert++; } if (bndVert != Gem_numBV(thee)) Vnm_print(2,"Gem_formChk: bndVert BAD\n"); /* go through all edges and check for consistency */ if (key >= 2) { Vnm_print(0,"Gem_formChk: ..EE..\n"); edgeOrderProb = 0; for (i=0; i<Gem_numEE(thee); i++) { eg = Gem_EE(thee,i); if ( (i>0) && (i % VPRTKEY) == 0 ) Vnm_print(0,"[CE:%d]",i); VJMPERR1( !Vnm_sigInt() ); /* check basic data */ /* (edge re-orderings ok; may not be error) */ if ((int)EE_id(eg) != i) edgeOrderProb = 1; if ( EE_vertex(eg,0) == VNULL ) Vnm_print(2,"Gem_formChk: vertex0 BAD edge <%d>\n",i); if ( EE_vertex(eg,1) == VNULL ) Vnm_print(2,"Gem_formChk: vertex1 BAD edge <%d>\n",i); /* check edge ring consistencies */ if ( ! VV_edgeInRing(EE_vertex(eg,0), eg) ) Vnm_print(2,"Gem_formChk: edgeInRing0 BAD edge <%d>\n",i); if ( ! VV_edgeInRing(EE_vertex(eg,1), eg) ) Vnm_print(2,"Gem_formChk: edgeInRing1 BAD edge <%d>\n",i); /* check midpoint data; should be VNULL */ if ( EE_midPoint(eg) != VNULL ) Vnm_print(2,"Gem_formChk: midPoint BAD edge <%d>\n",i); } if (edgeOrderProb) Vnm_print(0,"Gem_formChk: WARNING noncons edge order..)\n"); } /* go through all simplices and check for consistency */ Vnm_print(0,"Gem_formChk: ..SS..\n"); bndFace = 0; for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); if ( (i>0) && (i % VPRTKEY) == 0 ) Vnm_print(0,"[CS:%d]",i); VJMPERR1( !Vnm_sigInt() ); if ( (int)SS_id(sm) != i ) Vnm_print(2,"Gem_formChk: SS_id BAD sim <%d>\n",i); /* check for well-ordering of vertices through positive volume */ svol = Gem_simplexVolume(thee,sm); if (svol <= VSMALL) Vnm_print(2,"Gem_formChk: Degenerate sim <%d> has volume <%g>\n", i, svol); /* another check for well-ordering of vertices via positive volume */ if ( !Gem_orient(thee,sm) ) Vnm_print(2,"Gem_formChk: sim <%d> is badly ordered\n", i); /* no simplices should be marked for refinement/unrefinement */ if ( SS_refineKey(sm,0) != 0 ) Vnm_print(2,"Gem_formChk: SS_refineKey BAD sim <%d>\n",i); if ( SS_refineKey(sm,1) != 0 ) Vnm_print(2,"Gem_formChk: SS_refineKey BAD sim <%d>\n",i); /* get simplex face info */ for (j=0; j<4; j++) fType[j] = SS_faceType(sm,j); /* check vertex data and simplex ring structure */ for (j=0; j<Gem_dimVV(thee); j++) { v[j] = SS_vertex(sm,j); if ( v[j] == VNULL ) Vnm_print(2,"Gem_formChk: v[%d] BAD sim <%d>\n",j,i); if (key >= 1) { if ( ! VV_simplexInRing( v[j], sm ) ) Vnm_print(2,"Gem_formChk: sInR BAD sim <%d>\n",i); } } /* check simplex ring, boundary faces count, vertex type */ for (j=0; j<Gem_dimVV(thee); j++) { /* boundary vertex check */ if ( VBOUNDARY(fType[j]) ) { bndFace++; for (k=1; k<Gem_dimVV(thee); k++) { l=(j+k) % Gem_dimVV(thee); if ( VINTERIOR( VV_type(SS_vertex(sm,l)) ) ) Vnm_print(2,"Gem_formChk: Dv[%d] BAD sim <%d>\n", l,i); } } /* face check with all nabors; verify we are conforming */ if (key >= 3) { k=(j+1) % Gem_dimVV(thee); l=(k+1) % Gem_dimVV(thee); m=(l+1) % Gem_dimVV(thee); conformCount = 0; for (sm0=VV_firstSS(v[k]); sm0!=VNULL;sm0=SS_link(sm0,v[k])) { for (sm1=VV_firstSS(v[l]);sm1!=VNULL;sm1=SS_link(sm1,v[l])){ if (Gem_dim(thee) == 2) { if ((sm0!=sm) && (sm0==sm1)) conformCount++; } else { for (sm2=VV_firstSS(v[m]); sm2!=VNULL; sm2=SS_link(sm2,v[m])) { if ((sm0!=sm) && (sm0==sm1) && (sm0==sm2)) { conformCount++; } } } } } if ( VBOUNDARY(fType[j]) ) { if ( conformCount != 0 ) { Vnm_print(2,"Gem_formChk: conform BAD sim <%d>",i); Vnm_print(2,"..face <%d> should NOT have a nabor\n",j); } } else { if ( conformCount < 1 ) { Vnm_print(2,"Gem_formChk: conform BAD sim <%d>",i); Vnm_print(2,"..face <%d> should have a nabor\n",j); } if ( conformCount > 1 ) { Vnm_print(2,"Gem_formChk: conform BAD sim <%d>",i); Vnm_print(2,"..face <%d> has more than one nabor\n",j); } } } } } if ( bndFace != Gem_numBF(thee) ) Vnm_print(2,"Gem_formChk: numBF count BAD\n"); Vnm_print(0,"Gem_formChk: ..done.\n"); /* return with no errors */ Vnm_tstop(80, "form check"); return; VERROR1: Vnm_print(0,"Gem_formChk: ..done. [Early termination was forced]\n"); Vnm_tstop(80, "form check"); return; }