/* * *************************************************************************** * Routine: Aprx_internalBoundaries * * Purpose: Make internal boundary datastructures * * Author: Michael Holst and Daniel Reynolds * *************************************************************************** */ VPUBLIC int Aprx_internalBoundaries(Aprx *thee) { /* DAN-MIKE: */ int i, j, k, faceIDs, currID; SS *sm, *sm0; faceIDs = 0; for (i=0; i<Gem_numSS(thee->gm); i++) { sm = Gem_SS(thee->gm,i); for (j=0; j<Gem_dimVV(thee->gm); j++) { sm0 = SS_nabor(sm,j); if (SS_chart(sm0) != SS_chart(sm)) { for (k=0; k<Gem_dimVV(thee->gm); k++) { if (( sm0->d.fPtr[k] != VNULL ) && ( SS_nabor(sm0,k) == sm )) { currID = ((FF*)(sm0->d.fPtr[k]))->g.uid; } else { currID = faceIDs; faceIDs++; } } /* * allocate: sm->fptr->FF for face j of sm * FF_ctor needs: grandParentFace,sPtr,color */ sm->d.fPtr[j] = FF_ctor( currID, sm, SS_chart(sm0) ); } } } return 0; }
/* * *************************************************************************** * Routine: Gem_clearEdges * * Purpose: Clear all of the edge numbers in each simplex. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_clearEdges(Gem *thee) { # if defined(VG_ELEMENT) int i,j; SS *sm; # endif Vnm_tstart(70, "edge clear"); Vnm_print(0,"Gem_clearEdges: clearing edge numbers.."); /* go through simplices and clear edge numbers */ #if defined(VG_ELEMENT) for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); for (j=0; j<Gem_dimEE(thee); j++) { SS_setEdgeNumber(sm,j,-1); } } #endif Vnm_print(0,"..done.\n"); Vnm_tstop(70, "edge clear"); /* set the edge counter and return */ Gem_setNumVirtEE(thee, 0); }
/* /////////////////////////////////////////////////////////////////////////// // Routine: Vpee_numSS // // Author: Nathan Baker /////////////////////////////////////////////////////////////////////////// */ VPUBLIC int Vpee_numSS(Vpee *thee) { int num = 0; int isimp; for (isimp=0; isimp<Gem_numSS(thee->gm); isimp++) { if (SS_chart(Gem_SS(thee->gm, isimp)) == thee->localPartID) num++; } return num; }
/* * *************************************************************************** * Routine: Aprx_partSet * * Purpose: Clear the partitioning (set all colors to pcolor). * * Author: Michael Holst * *************************************************************************** */ VPUBLIC int Aprx_partSet(Aprx *thee, int pcolor) { int i; Vnm_print(0,"Aprx_partSet: [pc=%d] resetting partitioning..", pcolor); /* clear the mesh coloring */ for (i=0; i<Gem_numSS(thee->gm); i++) { SS_setChart(Gem_SS(thee->gm,i), pcolor); } Vnm_print(0,"..done.\n"); return 0; }
/* * *************************************************************************** * 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); }
/* * *************************************************************************** * Routine: Gem_countEdges * * Purpose: Count up all of the edges without actually creating them, and * keep track of the global edge numbers in each element. * * Notes: Based on a simplex traversal. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_countEdges(Gem *thee) { int edgeCnt; # if defined(VG_ELEMENT) int i,j,k,l,iv,iv2,edgeNum; SS *sm, *sm0, *sm2; VV *vx, *vx2; # endif Vnm_tstart(70, "edge count"); Vnm_print(0,"Gem_countEdges: counting edges (simplex traversal).."); /* go through simplices and count edges */ edgeCnt = 0; #if defined(VG_ELEMENT) for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); for (j=0; j<Gem_dimEE(thee); j++) { k = vmapEI[j][0]; l = vmapEI[j][1]; edgeNum = vmapE[k][l]; if (SS_edgeNumber(sm,edgeNum) < 0) { vx = SS_vertex(sm,k); vx2 = SS_vertex(sm,l); /* tell all simplices using edge what its number is */ for (sm2=VV_firstSS(vx2);sm2!=VNULL;sm2=SS_link(sm2,vx2)){ for (sm0=VV_firstSS(vx);sm0!=VNULL;sm0=SS_link(sm0,vx)){ if (sm0 == sm2) { iv = SS_vptr2localVnum(sm0, vx); iv2 = SS_vptr2localVnum(sm0, vx2); edgeNum = vmapE[iv][iv2]; SS_setEdgeNumber(sm0,edgeNum,edgeCnt); } } } edgeCnt++; } } } #endif Vnm_print(0,"..done. [edges=<%d>]\n",edgeCnt); Vnm_tstop(70, "edge count"); /* set the edge counter and return */ Gem_setNumVirtEE(thee, edgeCnt); }
/* * *************************************************************************** * Routine: Gem_createSimplexRings * * Purpose: Create all of the simplex rings. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_createSimplexRings(Gem *thee) { int i; SS *sm; Vnm_tstart(70, "simplex ring create"); Vnm_print(0,"Gem_createSimplexRings: creating simplex rings.."); /* go through all simplices and create simplex rings */ for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); SS_buildRing( sm ); } Vnm_print(0,"..done.\n"); Vnm_tstop(70, "simplex ring create"); }
/* * *************************************************************************** * Routine: Gem_createEdges * * Purpose: Create all of the edges. * * Notes: Based on a simplex traversal. * * We also set the edge numbers in the simplices while we we * are doing the edge creation. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_createEdges(Gem *thee) { int i,j,k,l,iDid,edgeCnt; # if defined(VG_ELEMENT) int edgeNum; EE *eg; # endif SS *sm; VV *v0, *v1; Vnm_tstart(70, "edge create"); Vnm_print(0,"Gem_createEdges: creating edges (simplex traversal).."); /* go through simplices and create edges if they don't already exist */ if (Gem_numEE(thee) == 0) { for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); for (j=0; j<Gem_dimEE(thee); j++) { k = vmapEI[j][0]; l = vmapEI[j][1]; v0 = SS_vertex(sm,k); v1 = SS_vertex(sm,l); # if defined(VG_ELEMENT) eg = Gem_findOrCreateEdge(thee, v0, v1, &iDid); edgeNum = EE_id( eg ); SS_setEdgeNumber(sm,j,edgeNum); # else (void*)Gem_findOrCreateEdge(thee, v0, v1, &iDid); # endif } } } edgeCnt = Gem_numEE(thee); Vnm_print(0,"..done. [edges=<%d>]\n",edgeCnt); Vnm_tstop(70, "edge create"); /* set the edge counter and return */ Gem_setNumVirtEE(thee, edgeCnt); }
/* * *************************************************************************** * Routine: MCsh_memChk * * Purpose: Print the exact current malloc usage. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void MCsh_memChk(MCsh *thee) { Vnm_print(0,"mc_vertex = %d;\n", Gem_numVV(thee->gm)); Vnm_print(0,"mc_simplex = %d;\n", Gem_numSS(thee->gm)); Vnm_print(0,"mc_memory = [\n"); Vnm_print(0,"%% --------------------------------------" "--------------------------------------\n"); Vnm_print(0,"%% Footprint Areas Malloc Free" " Highwater Class\n"), Vnm_print(0,"%% --------------------------------------" "--------------------------------------\n"); if (thee->gm != VNULL) Gem_memChk(thee->gm); if (thee->am != VNULL) AM_memChk(thee->am); Vmem_print(VNULL); Vnm_print(0,"%% --------------------------------------" "--------------------------------------\n"); Vmem_printTotal(); Vnm_print(0,"%% --------------------------------------" "--------------------------------------\n"); Vnm_print(0,"];\n"); if (thee->gm != VNULL) Gem_memChkMore(thee->gm); }
/* * *************************************************************************** * Routine: Gem_destroySimplexRings * * Purpose: Destroy all of the simplex rings. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_destroySimplexRings(Gem *thee) { int i; VV *vx; SS *sm; Vnm_print(0,"Gem_destroySimplexRings: destroying simplex rings.."); /* go through all vertices and reset firstSS pointer */ for (i=0; i<Gem_numVV(thee); i++) { vx = Gem_VV(thee,i); VV_setFirstSS( vx, VNULL ); } /* go through all simplices and reinitialize simplex rings */ for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); SS_initRing( sm ); } Vnm_print(0,"..done.\n"); }
/* * *************************************************************************** * Routine: Gem_countChk * * Purpose: Count all vertices, edges, faces, simplices, and do it * is cheaply as possible. Also do a form check. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_countChk(Gem *thee) { /* count vertices */ Gem_setNumVirtVV( thee, Gem_numVV(thee) ); /* count edges */ if ( Gem_numEE(thee) > 0 ) { /* edges DO exist; just use the current edge count */ Gem_setNumVirtEE( thee, Gem_numEE(thee) ); } else { /* edges DO NOT exist; create/destroy to get edge nums (fastest way) */ Gem_createEdges( thee ); Gem_destroyEdges( thee ); } /* count faces -- COSTLY! Only do if you have face DOF in the element */ Gem_countFaces( thee ); /* count simplices */ Gem_setNumVirtSS( thee, Gem_numSS(thee) ); /* check ringed vertex: 0=[v+s],1=[.+sRing],2=[.+eRing],3=[.+conform] */ /* Gem_formChk(thee,2); */ }
/* * *************************************************************************** * 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; }
VPUBLIC void Vcsm_init(Vcsm *thee) { /* Counters */ int iatom, jatom, isimp, jsimp, gotSimp; /* Atomic information */ Vatom *atom; double *position; /* Simplex/Vertex information */ SS *simplex; /* Basis function values */ if (thee == VNULL) { Vnm_print(2, "Vcsm_init: Error! Got NULL thee!\n"); VASSERT(0); } if (thee->gm == VNULL) { VASSERT(thee->gm != VNULL); Vnm_print(2, "Vcsm_init: Error! Got NULL thee->gm!\n"); VASSERT(0); } thee->nsimp = Gem_numSS(thee->gm); if (thee->nsimp <= 0) { Vnm_print(2, "Vcsm_init: Error! Got %d simplices!\n", thee->nsimp); VASSERT(0); } thee->natom = Valist_getNumberAtoms(thee->alist); /* Allocate and initialize space for the first dimensions of the * simplex-charge map, the simplex array, and the counters */ thee->sqm = Vmem_malloc(thee->vmem, thee->nsimp, sizeof(int *)); VASSERT(thee->sqm != VNULL); thee->nsqm = Vmem_malloc(thee->vmem, thee->nsimp, sizeof(int)); VASSERT(thee->nsqm != VNULL); for (isimp=0; isimp<thee->nsimp; isimp++) (thee->nsqm)[isimp] = 0; /* Count the number of charges per simplex. */ for (iatom=0; iatom<thee->natom; iatom++) { atom = Valist_getAtom(thee->alist, iatom); position = Vatom_getPosition(atom); gotSimp = 0; for (isimp=0; isimp<thee->nsimp; isimp++) { simplex = Gem_SS(thee->gm, isimp); if (Gem_pointInSimplex(thee->gm, simplex, position)) { (thee->nsqm)[isimp]++; gotSimp = 1; } } } /* Allocate the space for the simplex-charge map */ for (isimp=0; isimp<thee->nsimp; isimp++) { if ((thee->nsqm)[isimp] > 0) { thee->sqm[isimp] = Vmem_malloc(thee->vmem, (thee->nsqm)[isimp], sizeof(int)); VASSERT(thee->sqm[isimp] != VNULL); } } /* Finally, set up the map */ for (isimp=0; isimp<thee->nsimp; isimp++) { jsimp = 0; simplex = Gem_SS(thee->gm, isimp); for (iatom=0; iatom<thee->natom; iatom++) { atom = Valist_getAtom(thee->alist, iatom); position = Vatom_getPosition(atom); /* Check to see if the atom's in this simplex */ if (Gem_pointInSimplex(thee->gm, simplex, position)) { /* Assign the entries in the next vacant spot */ (thee->sqm)[isimp][jsimp] = iatom; jsimp++; } } } thee->msimp = thee->nsimp; /* Allocate space for the charge-simplex map */ thee->qsm = Vmem_malloc(thee->vmem, thee->natom, sizeof(int *)); VASSERT(thee->qsm != VNULL); thee->nqsm = Vmem_malloc(thee->vmem, thee->natom, sizeof(int)); VASSERT(thee->nqsm != VNULL); for (iatom=0; iatom<thee->natom; iatom++) (thee->nqsm)[iatom] = 0; /* Loop through the list of simplices and count the number of times * each atom appears */ for (isimp=0; isimp<thee->nsimp; isimp++) { for (iatom=0; iatom<thee->nsqm[isimp]; iatom++) { jatom = thee->sqm[isimp][iatom]; thee->nqsm[jatom]++; } } /* Do a TIME-CONSUMING SANITY CHECK to make sure that each atom was * placed in at simplex */ for (iatom=0; iatom<thee->natom; iatom++) { if (thee->nqsm[iatom] == 0) { Vnm_print(2, "Vcsm_init: Atom %d not placed in simplex!\n", iatom); VASSERT(0); } } /* Allocate the appropriate amount of space for each entry in the * charge-simplex map and clear the counter for re-use in assignment */ for (iatom=0; iatom<thee->natom; iatom++) { thee->qsm[iatom] = Vmem_malloc(thee->vmem, (thee->nqsm)[iatom], sizeof(int)); VASSERT(thee->qsm[iatom] != VNULL); thee->nqsm[iatom] = 0; } /* Assign the simplices to atoms */ for (isimp=0; isimp<thee->nsimp; isimp++) { for (iatom=0; iatom<thee->nsqm[isimp]; iatom++) { jatom = thee->sqm[isimp][iatom]; thee->qsm[jatom][thee->nqsm[jatom]] = isimp; thee->nqsm[jatom]++; } } thee->initFlag = 1; }
/* * *************************************************************************** * Routine: Gem_formFix * * Purpose: Make some specified hacked fix to a given mesh. * * Notes: key==0 --> ? * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_formFix(Gem *thee, int key) { int i, j, k, l, m, nabors, btype; double radk, radl, radm, myTol; VV *v[4]; SS *sm, *sm0, *sm1, *sm2; /* input check and some i/o */ btype = key; VASSERT( (0 <= btype) && (btype <= 2) ); /* go through all simplices and zero all boundary faces */ Vnm_print(0,"Gem_makeBnd: zeroing boundary faces/vertices.."); Gem_setNumBF(thee, 0); Gem_setNumBV(thee, 0); for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); if ( (i>0) && (i % VPRTKEY) == 0 ) Vnm_print(0,"[BS:%d]",i); /* get local vertices */ for (j=0; j<Gem_dimVV(thee); j++) v[j] = SS_vertex(sm,j); /* reset all vertices and faces to interior type */ for (j=0; j<Gem_dimVV(thee); j++) { /* the other three local vertex/face numbers besides "j" */ k=(j+1) % Gem_dimVV(thee); l=(k+1) % Gem_dimVV(thee); m=(l+1) % Gem_dimVV(thee); SS_setFaceType(sm, j, 0); VV_setType(v[k], 0); VV_setType(v[l], 0); if (Gem_dim(thee) == 3) VV_setType(v[m], 0); } } Vnm_print(0,"..done.\n"); /* are we done */ /* if (btype == 0) return; */ /* okay now make a boundary */ Vnm_print(0,"Gem_makeBnd: rebuilding boundary faces/vertices.."); for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); if ( (i>0) && (i % VPRTKEY) == 0 ) Vnm_print(0,"[BS:%d]",i); /* get local vertices */ for (j=0; j<Gem_dimVV(thee); j++) v[j] = SS_vertex(sm,j); /* rebuild everything */ for (j=0; j<Gem_dimVV(thee); j++) { /* the other three local vertex/face numbers besides "j" */ k=(j+1) % Gem_dimVV(thee); l=(k+1) % Gem_dimVV(thee); m=(l+1) % Gem_dimVV(thee); /* look for a face nabor sharing face "j" (opposite vertex "j") */ nabors = 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)) nabors++; } else { for (sm2=VV_firstSS(v[m]); sm2!=VNULL; sm2=SS_link(sm2,v[m])) { if ((sm0!=sm) && (sm0==sm1) && (sm0==sm2)) { nabors++; } } } } } /* if no one there, then face "j" is actually a boundary face */ if (nabors == 0) { myTol = 1.0e-2; if ( ( VABS(VV_coord(v[k],2) - 0.0) < myTol) && ( VABS(VV_coord(v[l],2) - 0.0) < myTol) && ( VABS(VV_coord(v[m],2) - 0.0) < myTol) ) { btype = 1; } else if ( ( VABS(VV_coord(v[k],2) - 68.03512) < myTol) && ( VABS(VV_coord(v[l],2) - 68.03512) < myTol) && ( VABS(VV_coord(v[m],2) - 68.03512) < myTol) ) { btype = 3; } else { radk = VSQRT( VSQR( VV_coord(v[k],0) ) + VSQR( VV_coord(v[k],1) ) ); radl = VSQRT( VSQR( VV_coord(v[l],0) ) + VSQR( VV_coord(v[l],1) ) ); radm = VSQRT( VSQR( VV_coord(v[m],0) ) + VSQR( VV_coord(v[m],1) ) ); if ( ( VABS(radk - 1.5) < myTol) && ( VABS(radl - 1.5) < myTol) && ( VABS(radm - 1.5) < myTol) ) { btype = 2; } else if ( ( VABS(radk - 2.0) < myTol) && ( VABS(radl - 2.0) < myTol) && ( VABS(radm - 2.0) < myTol) ){ btype = 4; } else { btype = 0; } } SS_setFaceType(sm, j, btype); Gem_numBFpp(thee); if (VINTERIOR( VV_type(v[k])) ) { VV_setType(v[k], btype); Gem_numBVpp(thee); } if (VINTERIOR( VV_type(v[l])) ) { VV_setType(v[l], btype); Gem_numBVpp(thee); } if (Gem_dim(thee) == 3) { if (VINTERIOR( VV_type(v[m])) ) { VV_setType(v[m], btype); Gem_numBVpp(thee); } } } } } Vnm_print(0,"..done.\n"); }
/* * *************************************************************************** * Routine: Gem_makeBndExt * * Purpose: Mark selected boundary faces in a special way. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_makeBndExt(Gem *thee, int key) { int i, j, k, l, m, p, q, nabors, btype, done, btypeGeneric; VV *v[4]; SS *sm, *sm0, *sm1, *sm2; double x[4][3], xchk; /* go through all simplices and zero all boundary faces */ Vnm_print(0,"Gem_makeBnd: zeroing boundary faces/vertices.."); Gem_setNumBF(thee, 0); Gem_setNumBV(thee, 0); for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); if ( (i>0) && (i % VPRTKEY) == 0 ) Vnm_print(0,"[BS:%d]",i); /* get local vertices */ for (j=0; j<Gem_dimVV(thee); j++) v[j] = SS_vertex(sm,j); /* reset all vertices and faces to interior type */ for (j=0; j<Gem_dimVV(thee); j++) { /* the other three local vertex/face numbers besides "j" */ k=(j+1) % Gem_dimVV(thee); l=(k+1) % Gem_dimVV(thee); m=(l+1) % Gem_dimVV(thee); SS_setFaceType(sm, j, 0); VV_setType(v[k], 0); VV_setType(v[l], 0); if (Gem_dim(thee) == 3) VV_setType(v[m], 0); } } Vnm_print(0,"..done.\n"); /* okay now make a boundary */ Vnm_print(0,"Gem_makeBnd: rebuilding boundary faces/vertices.."); for (i=0; i<Gem_numSS(thee); i++) { sm = Gem_SS(thee,i); if ( (i>0) && (i % VPRTKEY) == 0 ) Vnm_print(0,"[BS:%d]",i); /* get local vertices */ for (j=0; j<Gem_dimVV(thee); j++) v[j] = SS_vertex(sm,j); /* rebuild everything */ for (j=0; j<Gem_dimVV(thee); j++) { /* the other three local vertex/face numbers besides "j" */ k=(j+1) % Gem_dimVV(thee); l=(k+1) % Gem_dimVV(thee); m=(l+1) % Gem_dimVV(thee); /* look for a face nabor sharing face "j" (opposite vertex "j") */ nabors = 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)) nabors++; } else { for (sm2=VV_firstSS(v[m]); sm2!=VNULL; sm2=SS_link(sm2,v[m])) { if ((sm0!=sm) && (sm0==sm1) && (sm0==sm2)) { nabors++; } } } } } /* if no one there, then face "j" is actually a boundary face */ if (nabors == 0) { /* grab coordinates of the vertices of this face */ for (q=0; q<Gem_dim(thee); q++) { x[0][q] = VV_coord(v[k],q); } for (q=0; q<Gem_dim(thee); q++) { x[1][q] = VV_coord(v[l],q); } if (Gem_dim(thee) == 3) { for (q=0; q<Gem_dim(thee); q++) { x[2][q] = VV_coord(v[m],q); } } /* default is interior; should not occur! */ btypeGeneric = 18; done = 0; btype = btypeGeneric; /* ---------- check for base marking ---------- */ xchk = 0.0; for (p=0; p<Gem_dim(thee); p++) { xchk += VABS( x[p][1] - (-1.0) ); } if (xchk < VSMALL) { done = 1; btype = 1; } /* ---------- check for base marking again ---------- */ xchk = 0.0; for (p=0; p<Gem_dim(thee); p++) { xchk += VABS( x[p][1] - ( 0.0) ); } if (xchk < VSMALL) { done = 1; btype = 18; } /* ---------- check for first section ---------- */ if (!done) { done = 1; btype = 2; for (p=0; p<Gem_dim(thee); p++) { if (! ( ( 1.9 <= x[p][0]) && ( 6.1 >= x[p][0]) && (-VSMALL <= x[p][1]) && (-1.1 <= x[p][2]) && ( 1.1 >= x[p][2]) )) { done = 0; btype = btypeGeneric; } } if (done) { xchk = 0.0; for (p=0; p<Gem_dim(thee); p++) { xchk += VABS( x[p][1] - 10.0 ); } if (xchk < VSMALL) { btype = 10; } } } /* ---------- check for second section ---------- */ if (!done) { done = 1; btype = 4; for (p=0; p<Gem_dim(thee); p++) { if (! ( ( 7.9 <= x[p][0]) && (12.1 >= x[p][0]) && (-VSMALL <= x[p][1]) && (-1.1 <= x[p][2]) && ( 1.1 >= x[p][2]) )) { done = 0; btype = btypeGeneric; } } if (done) { xchk = 0.0; for (p=0; p<Gem_dim(thee); p++) { xchk += VABS( x[p][1] - 10.0 ); } if (xchk < VSMALL) { btype = 12; } } } /* ---------- check for third section ---------- */ if (!done) { done = 1; btype = 6; for (p=0; p<Gem_dim(thee); p++) { if (! ( (13.9 <= x[p][0]) && (18.1 >= x[p][0]) && (-VSMALL <= x[p][1]) && (-1.1 <= x[p][2]) && ( 1.1 >= x[p][2]) )) { done = 0; btype = btypeGeneric; } } if (done) { xchk = 0.0; for (p=0; p<Gem_dim(thee); p++) { xchk += VABS( x[p][1] - 10.0 ); } if (xchk < VSMALL) { btype = 14; } } } /* ---------- check for fourth section ---------- */ if (!done) { done = 1; btype = 8; for (p=0; p<Gem_dim(thee); p++) { if (! ( (19.9 <= x[p][0]) && (24.1 >= x[p][0]) && (-VSMALL <= x[p][1]) && (-1.1 <= x[p][2]) && ( 1.1 >= x[p][2]) )) { done = 0; btype = btypeGeneric; } } if (done) { xchk = 0.0; for (p=0; p<Gem_dim(thee); p++) { xchk += VABS( x[p][1] - 10.0 ); } if (xchk < VSMALL) { btype = 16; } } } /* should have been marked with SOME boundary type */ VASSERT( 0 != btype ); /* set the facetype */ SS_setFaceType(sm, j, btype); Gem_numBFpp(thee); /* set the vertex types (dirichlet overrides robin) */ if (!VDIRICHLET( VV_type(v[k])) ) { if (VINTERIOR( VV_type(v[k])) ) { Gem_numBVpp(thee); } VV_setType(v[k], btype); } if (!VDIRICHLET( VV_type(v[l])) ) { if (VINTERIOR( VV_type(v[l])) ) { Gem_numBVpp(thee); } VV_setType(v[l], btype); } if (Gem_dim(thee) == 3) { if (!VDIRICHLET( VV_type(v[m])) ) { if (VINTERIOR( VV_type(v[m])) ) { Gem_numBVpp(thee); } VV_setType(v[m], btype); } } } } } Vnm_print(0,"..done.\n"); }
/* /////////////////////////////////////////////////////////////////////////// // Routine: Vpee_markRefine // // Author: Nathan Baker (and Michael Holst: the author of AM_markRefine, on // which this is based) /////////////////////////////////////////////////////////////////////////// */ VPUBLIC int Vpee_markRefine(Vpee *thee, AM *am, int level, int akey, int rcol, double etol, int bkey ) { Aprx *aprx; int marked = 0, markMe, i, smid, count, currentQ; double minError = 0.0, maxError = 0.0, errEst = 0.0, mlevel, barrier; SS *sm; VASSERT(thee != VNULL); /* Get the Aprx object from AM */ aprx = am->aprx; /* input check and some i/o */ if ( ! ((-1 <= akey) && (akey <= 4)) ) { Vnm_print(0,"Vpee_markRefine: bad refine key; simplices marked = %d\n", marked); return marked; } /* For uniform markings, we have no effect */ if ((-1 <= akey) && (akey <= 0)) { marked = Gem_markRefine(thee->gm, akey, rcol); return marked; } /* Informative I/O */ if (akey == 2) { Vnm_print(0,"Vpee_estRefine: using Aprx_estNonlinResid().\n"); } else if (akey == 3) { Vnm_print(0,"Vpee_estRefine: using Aprx_estLocalProblem().\n"); } else if (akey == 4) { Vnm_print(0,"Vpee_estRefine: using Aprx_estDualProblem().\n"); } else { Vnm_print(0,"Vpee_estRefine: bad key given; simplices marked = %d\n", marked); return marked; } if (thee->killFlag == 0) { Vnm_print(0, "Vpee_markRefine: No error attenuation -- simplices in all partitions will be marked.\n"); } else if (thee->killFlag == 1) { Vnm_print(0, "Vpee_markRefine: Maximum error attenuation -- only simplices in local partition will be marked.\n"); } else if (thee->killFlag == 2) { Vnm_print(0, "Vpee_markRefine: Spherical error attenutation -- simplices within a sphere of %4.3f times the size of the partition will be marked\n", thee->killParam); } else if (thee->killFlag == 2) { Vnm_print(0, "Vpee_markRefine: Neighbor-based error attenuation -- simplices in the local and neighboring partitions will be marked [NOT IMPLEMENTED]!\n"); VASSERT(0); } else { Vnm_print(2,"Vpee_markRefine: bogus killFlag given; simplices marked = %d\n", marked); return marked; } /* set the barrier type */ mlevel = (etol*etol) / Gem_numSS(thee->gm); if (bkey == 0) { barrier = (etol*etol); Vnm_print(0,"Vpee_estRefine: forcing [err per S] < [TOL] = %g\n", barrier); } else if (bkey == 1) { barrier = mlevel; Vnm_print(0,"Vpee_estRefine: forcing [err per S] < [(TOL^2/numS)^{1/2}] = %g\n", VSQRT(barrier)); } else { Vnm_print(0,"Vpee_estRefine: bad bkey given; simplices marked = %d\n", marked); return marked; } /* timer */ Vnm_tstart(30, "error estimation"); /* count = num generations to produce from marked simplices (minimally) */ count = 1; /* must be >= 1 */ /* check the refinement Q for emptyness */ currentQ = 0; if (Gem_numSQ(thee->gm,currentQ) > 0) { Vnm_print(0,"Vpee_markRefine: non-empty refinement Q%d....clearing..", currentQ); Gem_resetSQ(thee->gm,currentQ); Vnm_print(0,"..done.\n"); } if (Gem_numSQ(thee->gm,!currentQ) > 0) { Vnm_print(0,"Vpee_markRefine: non-empty refinement Q%d....clearing..", !currentQ); Gem_resetSQ(thee->gm,!currentQ); Vnm_print(0,"..done.\n"); } VASSERT( Gem_numSQ(thee->gm,currentQ) == 0 ); VASSERT( Gem_numSQ(thee->gm,!currentQ) == 0 ); /* clear everyone's refinement flags */ Vnm_print(0,"Vpee_markRefine: clearing all simplex refinement flags.."); for (i=0; i<Gem_numSS(thee->gm); i++) { if ( (i>0) && (i % VPRTKEY) == 0 ) Vnm_print(0,"[MS:%d]",i); sm = Gem_SS(thee->gm,i); SS_setRefineKey(sm,currentQ,0); SS_setRefineKey(sm,!currentQ,0); SS_setRefinementCount(sm,0); } Vnm_print(0,"..done.\n"); /* NON-ERROR-BASED METHODS */ /* Simplex flag clearing */ if (akey == -1) return marked; /* Uniform & user-defined refinement*/ if ((akey == 0) || (akey == 1)) { smid = 0; while ( smid < Gem_numSS(thee->gm)) { /* Get the simplex and find out if it's markable */ sm = Gem_SS(thee->gm,smid); markMe = Vpee_ourSimp(thee, sm, rcol); if (markMe) { if (akey == 0) { marked++; Gem_appendSQ(thee->gm,currentQ, sm); SS_setRefineKey(sm,currentQ,1); SS_setRefinementCount(sm,count); } else if (Vpee_userDefined(thee, sm)) { marked++; Gem_appendSQ(thee->gm,currentQ, sm); SS_setRefineKey(sm,currentQ,1); SS_setRefinementCount(sm,count); } } smid++; } } /* ERROR-BASED METHODS */ /* gerror = global error accumulation */ aprx->gerror = 0.; /* traverse the simplices and process the error estimates */ Vnm_print(0,"Vpee_markRefine: estimating error.."); smid = 0; while ( smid < Gem_numSS(thee->gm)) { /* Get the simplex and find out if it's markable */ sm = Gem_SS(thee->gm,smid); markMe = Vpee_ourSimp(thee, sm, rcol); if ( (smid>0) && (smid % VPRTKEY) == 0 ) Vnm_print(0,"[MS:%d]",smid); /* Produce an error estimate for this element if it is in the set */ if (markMe) { if (akey == 2) { errEst = Aprx_estNonlinResid(aprx, sm, am->u,am->ud,am->f); } else if (akey == 3) { errEst = Aprx_estLocalProblem(aprx, sm, am->u,am->ud,am->f); } else if (akey == 4) { errEst = Aprx_estDualProblem(aprx, sm, am->u,am->ud,am->f); } VASSERT( errEst >= 0. ); /* if error estimate above tol, mark element for refinement */ if ( errEst > barrier ) { marked++; Gem_appendSQ(thee->gm,currentQ, sm); /*add to refinement Q*/ SS_setRefineKey(sm,currentQ,1); /* note now on refine Q */ SS_setRefinementCount(sm,count); /* refine X many times? */ } /* keep track of min/max errors over the mesh */ minError = VMIN2( VSQRT(VABS(errEst)), minError ); maxError = VMAX2( VSQRT(VABS(errEst)), maxError ); /* store the estimate */ Bvec_set( aprx->wev, smid, errEst ); /* accumlate into global error (errEst is SQUAREd already) */ aprx->gerror += errEst; /* otherwise store a zero for the estimate */ } else { Bvec_set( aprx->wev, smid, 0. ); } smid++; } /* do some i/o */ Vnm_print(0,"..done. [marked=<%d/%d>]\n",marked,Gem_numSS(thee->gm)); Vnm_print(0,"Vpee_estRefine: TOL=<%g> Global_Error=<%g>\n", etol, aprx->gerror); Vnm_print(0,"Vpee_estRefine: (TOL^2/numS)^{1/2}=<%g> Max_Ele_Error=<%g>\n", VSQRT(mlevel),maxError); Vnm_tstop(30, "error estimation"); /* check for making the error tolerance */ if ((bkey == 1) && (aprx->gerror <= etol)) { Vnm_print(0, "Vpee_estRefine: *********************************************\n"); Vnm_print(0, "Vpee_estRefine: Global Error criterion met; setting marked=0.\n"); Vnm_print(0, "Vpee_estRefine: *********************************************\n"); marked = 0; } /* return */ return marked; }
/* * *************************************************************************** * Routine: Aprx_partOne * * Purpose: Do a single bisection step. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC int Aprx_partOne(Aprx *thee, int pkey, int pwht, int pcolor, int poff) { int i, j, k, ii, dim, rc, numC, numS, sCount1, sCount2, cutOff; int *pord, *ford, *rord; double tm, cm[3], *evec, eCount1, eCount2, ecutOff, ePart; simHelper *simH; SS *sm; Vnm_print(0,"Aprx_partOne: [pcolor=%d] partitioning:\n", pcolor); /* dimensions */ dim = Gem_dim(thee->gm); numS = Gem_numSS(thee->gm); /* how many guys of this color */ numC = 0; for (i=0; i<numS; i++) { sm = Gem_SS(thee->gm,i); if ((int)SS_chart(sm) == pcolor) numC++; } /* grab some temporary storage */ pord = Vmem_malloc( thee->vmem, numC, sizeof(int) ); ford = Vmem_malloc( thee->vmem, numC, sizeof(int) ); rord = Vmem_malloc( thee->vmem, numS, sizeof(int) ); evec = Vmem_malloc( thee->vmem, numC, sizeof(double) ); simH = Vmem_malloc( thee->vmem, numC, sizeof(simHelper) ); /* build the simplex helper vector and forward/reverse ordering arrays */ ePart = 0.; tm = 0.; for (i=0; i<3; i++) cm[i] = 0.; i = 0; for (ii=0; ii<Gem_numSS(thee->gm); ii++) { sm = Gem_SS(thee->gm,ii); if ((int)SS_chart(sm) == pcolor) { /* initialization */ pord[i] = i; ford[i] = ii; rord[ii] = i; evec[i] = 0.; simH[i].color = SS_chart(sm); simH[i].diag = 0; simH[i].mass = 1.; simH[i].error = Bvec_valB( thee->wev, 0, ii ); for (j=0; j<4; j++) { simH[i].faceId[j] = -1; } /* accumulate error into total for this partition */ ePart += simH[i].error; /* baricenter of this simplex and center of mass of all */ tm += simH[i].mass; for (j=0; j<3; j++) { simH[i].bc[j] = 0.; for (k=0; k<dim+1; k++) { simH[i].bc[j] += VV_coord(SS_vertex(sm,k),j); } simH[i].bc[j] /= (dim+1.); cm[j] += (simH[i].mass * simH[i].bc[j]); } i++; } else { rord[ii] = -1; } } VASSERT( i == numC ); for (i=0; i<3; i++) cm[i] /= tm; /* translate coordinate systems so center of mass is at origin */ for (i=0; i<numC; i++) { for (j=0; j<3; j++) { simH[i].bc[j] -= cm[j]; } } /* okay partition it */ if (pkey == 0) { rc = Aprx_partInert(thee, pcolor, numC, evec, simH); } else if (pkey == 1) { rc = Aprx_partSpect(thee, pcolor, numC, evec, simH, ford, rord, 0); } else if (pkey == 2) { rc = Aprx_partSpect(thee, pcolor, numC, evec, simH, ford, rord, 1); } else if (pkey == 3) { rc = Aprx_partInert(thee, pcolor, numC, evec, simH); rc = Aprx_partSpect(thee, pcolor, numC, evec, simH, ford, rord, 0); } else { Vnm_print(2,"Aprx_partOne: illegal pkey value of <%d>\n", pkey); rc = -1; } /* sort the values small-to-big; allows us to have equi-size submeshes */ Vnm_dqsortOrd(evec, pord, numC); /* color the mesh using the ordering given by the eigenvector */ sCount1 = 0; sCount2 = 0; eCount1 = 0.; eCount2 = 0.; if (rc >= 0) { /* weighted partitioning; divides into two sets of equal error */ if (pwht == 1) { ecutOff = ( ePart / 2. ); for (i=0; i<numC; i++) { sm = Gem_SS(thee->gm,ford[pord[i]]); if (eCount1 < ecutOff) { SS_setChart( sm, pcolor ); sCount1++; eCount1 += Bvec_valB( thee->wev, 0, ford[pord[i]] ); } else { SS_setChart( sm, pcolor+poff ); sCount2++; eCount2 += Bvec_valB( thee->wev, 0, ford[pord[i]] ); } } /* unweighted partitioning; divides into two sets of equal number */ } else { cutOff = numC / 2; for (i=0; i<numC; i++) { sm = Gem_SS(thee->gm,ford[pord[i]]); if (sCount1 < cutOff) { SS_setChart( sm, pcolor ); sCount1++; eCount1 += Bvec_valB( thee->wev, 0, ford[pord[i]] ); } else { SS_setChart( sm, pcolor+poff ); sCount2++; eCount2 += Bvec_valB( thee->wev, 0, ford[pord[i]] ); } } } } /* free the temporary storage */ Vmem_free( thee->vmem, numC, sizeof(int), (void**)&pord ); Vmem_free( thee->vmem, numC, sizeof(int), (void**)&ford ); Vmem_free( thee->vmem, numS, sizeof(int), (void**)&rord ); Vmem_free( thee->vmem, numC, sizeof(double), (void**)&evec ); Vmem_free( thee->vmem, numC, sizeof(simHelper), (void**)&simH ); Vnm_print(0,"Aprx_partOne: done."); Vnm_print(0," [c1=%d,c2=%d,e1=%8.2e,e2=%8.2e,eT=%8.2e]\n", sCount1, sCount2, eCount1, eCount2, (eCount1+eCount2)); if (sCount1 == 0) { Vnm_print(2,"Aprx_partOne: ERROR: first partition has NO SIMPLICES!\n"); } if (sCount2 == 0) { Vnm_print(2,"Aprx_partOne: ERROR: second partition has NO SIMPLICES!\n"); } return rc; }
/* * *************************************************************************** * Routine: Aprx_partSmooth * * Purpose: Smooth the partitioning. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC int Aprx_partSmooth(Aprx *thee) { int i, j, dim, max, maxc, tcount, myc, mycount; int c[4], count[4]; SS *sm, *sm0; Vnm_print(0,"Aprx_partSmooth: smoothing partitioning.."); dim = Gem_dim(thee->gm); /* smooth the mesh coloring */ tcount = 0; for (i=0; i<Gem_numSS(thee->gm); i++) { sm = Gem_SS(thee->gm, i); myc = SS_chart(sm); mycount = 0; for (j=0; j<4; j++) { count[j] = 0; c[j] = -1; } for (j=0; j<Gem_dimVV(thee->gm); j++) { sm0 = SS_nabor(sm,j); if (sm0 != VNULL) { if ((int)SS_chart(sm0) == myc) { mycount++; } else { count[j] = 1; c[j] = SS_chart(sm0); } } } if ((c[0] >= 0) && (c[0] == c[1])) count[0]++; if ((c[0] >= 0) && (c[0] == c[2])) count[0]++; if ((c[0] >= 0) && (c[0] == c[3])) count[0]++; if ((c[1] >= 0) && (c[1] == c[0])) count[1]++; if ((c[1] >= 0) && (c[1] == c[2])) count[1]++; if ((c[1] >= 0) && (c[1] == c[3])) count[1]++; if ((c[2] >= 0) && (c[2] == c[0])) count[2]++; if ((c[2] >= 0) && (c[2] == c[1])) count[2]++; if ((c[2] >= 0) && (c[2] == c[3])) count[2]++; if (dim == 3) { if ((c[3] >= 0) && (c[3] == c[0])) count[3]++; if ((c[3] >= 0) && (c[3] == c[1])) count[3]++; if ((c[3] >= 0) && (c[3] == c[2])) count[3]++; } max = mycount; maxc = -1; for (j=0; j<Gem_dimVV(thee->gm); j++) { if (count[j] > max) { max = count[j]; maxc = j; } } if (max > mycount) { SS_setChart(sm, c[maxc]); tcount++; } } Vnm_print(0,"..done. [%d simplices changed color]\n", tcount); return tcount; }
/* * *************************************************************************** * 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: Gem_createAndInitSS * * Purpose: Create and initialize a new simplex; return a point to it. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC SS* Gem_createAndInitSS(Gem *thee) { SS *sm = Gem_createSS(thee); SS_init(sm, Gem_dim(thee), Gem_numSS(thee)-1 ); return sm; }