void apbsfinal_() { unsigned long int bytesTotal, highWater; int i; VASSERT(nosh != VNULL); /* Kill the saved potential Vgrids */ for (i=0; i<2; i++){ Vgrid_dtor(&permU[i]); Vgrid_dtor(&indU[i]); Vgrid_dtor(&nlIndU[i]); } Valist_dtor(&alist[0]); /* Saving the kappa and dielectric maps is under consideration. killKappaMaps(nosh, kappaMap); killDielMaps(nosh, dielXMap, dielYMap, dielZMap); */ NOsh_dtor(&nosh); /* Clean up MALOC structures */ bytesTotal = Vmem_bytesTotal(); highWater = Vmem_highWaterTotal(); /* printf(" Final APBS memory usage: %4.3f MB total, %4.3f MB high water\n\n", (double)(bytesTotal)/(1024.*1024.),(double)(highWater)/(1024.*1024.)); */ Vmem_dtor(&mem); Vnm_tstop(APBS_TIMER_WALL_CLOCK, "APBS WALL CLOCK"); Vcom_finalize(); Vcom_dtor(&com); }
/* * *************************************************************************** * 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: Bmat_galerkin * * Purpose: Enforce the Galerkin conditions algebraically. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Bmat_galerkin(Bmat *thee, Bmat *rmat, Bmat *amat, Bmat *pmat) { int p, q, numR, numC; char bname[15]; #if 0 Vnm_tstart(30, "triple matrix product"); Vnm_print(0,"Bmat_galerkin: triple matrix product.."); #endif /* go through blocks and form B=r*A*p */ for (p=0; p<thee->numB; p++) { for (q=0; q<thee->numB; q++) { /* form the product: B=R*A*P */ if ( !Bmat_mirror(thee,p,q) ) { /* hold onto row/col sizes */ numR = Bmat_numR(thee,p,q); numC = Bmat_numR(thee,p,q); /* reinitialize the block */ Mat_dtor(&thee->AD[p][q]); sprintf(bname, "%s_%d%d", thee->name, p, q); thee->AD[p][q] = Mat_ctor(thee->vmem, bname, numR, numC); /* don't forget the possible mirror!!! */ if ( Bmat_mirror(thee,q,p) ) { thee->AD[q][p] = thee->AD[p][q]; } /* now form the triple matrix product */ Mat_galerkin(thee->AD[p][q], rmat->AD[p][p], amat->AD[p][q], pmat->AD[q][q]); } } } /* * Now, insert dirichlet identity equations due to singular P. * * (The new matrix is singular otherwise, since the prolongation matrix * has zero COLUMNS at coarse dirichlet points and zero ROWS at fine * dirichlet points, with the reverse situation for the restriction.) */ Bmat_diri( thee ); /* Set the coarse and fine pointers within the matrices. */ thee->fine = amat; amat->coarse = thee; #if 0 Vnm_print(0,"..done.\n"); Vnm_tstop(30, "triple matrix product"); #endif }
/* * *************************************************************************** * 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: 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: 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; }
/* * *************************************************************************** * Routine: Gem_speedChk * * Purpose: Calculate the cost to traverse the various structures. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void Gem_speedChk(Gem *thee) { const int N = 100000; const int loop = 1000; int i, j, itotal; double rtotal; double *d0; int *i0; VV *v0, *vx; Vnm_print(0,"Gem_speedChk: Creating normal double vec...\n"); Vnm_tstart(90, "double vec creation"); d0 = Vmem_malloc( thee->vmem, N, sizeof(double) ); Vnm_tstop(90, "double vec creation"); Vnm_print(0,"Gem_speedChk: Creating normal int vec...\n"); Vnm_tstart(90, "int vec creation"); i0 = Vmem_malloc( thee->vmem, N, sizeof(int) ); Vnm_tstop(90, "int vec creation"); Vnm_print(0,"Gem_speedChk: Creating normal VV vec...\n"); Vnm_tstart(90, "VV vec creation"); v0 = Vmem_malloc( thee->vmem, N, sizeof(VV) ); Vnm_tstop(90, "VV vec creation"); Vnm_print(0,"Gem_speedChk: Creating Vset VV vec...\n"); Vnm_tstart(90, "VV Vset creation"); for (i=0; i<N; i++) vx = Gem_createAndInitVV(thee); Vnm_tstop(90, "VV Vset creation"); Vnm_print(0,"Gem_speedChk: Looping over double vec...\n"); Vnm_print(0,"Gem_speedChk: N, loop = %d, %d\n", N, loop); Vnm_tstart(90, "double vec loop"); rtotal = 0.0; for (j=0; j<loop; j++) for (i=0; i<N; i++) rtotal += ( j * d0[i] * (j+d0[i]) ); Vnm_tstop(90, "double vec loop"); Vnm_print(0,"Gem_speedChk: Looping over int vec...\n"); Vnm_print(0,"Gem_speedChk: N, loop = %d, %d\n", N, loop); Vnm_tstart(90, "int vec loop"); itotal = 0; for (j=0; j<loop; j++) for (i=0; i<N; i++) itotal += ( j * i0[i] * (j+i0[i]) ); Vnm_tstop(90, "int vec loop"); Vnm_print(0,"Gem_speedChk: Looping over VV vec...\n"); Vnm_print(0,"Gem_speedChk: N, loop = %d, %d\n", N, loop); Vnm_tstart(90, "VV vec loop"); itotal = 0; for (j=0; j<loop; j++) for (i=0; i<N; i++) itotal += ( j * VV_id(&v0[i]) * (j+VV_id(&v0[i])) ); Vnm_tstop(90, "VV vec loop"); Vnm_print(0,"Gem_speedChk: Looping (index) over Vset VV vec...\n"); Vnm_print(0,"Gem_speedChk: N, loop = %d, %d\n", N, loop); Vnm_tstart(90, "VV Vset loop"); itotal = 0; for (j=0; j<loop; j++) for (i=0; i<N; i++) itotal += ( j * VV_id(Gem_VV(thee,i)) * (j + VV_id(Gem_VV(thee,i)) ) ); Vnm_tstop(90, "VV Vset loop"); Vnm_print(0,"Gem_speedChk: Looping (iterator) over Vset VV vec...\n"); Vnm_print(0,"Gem_speedChk: N, loop = %d, %d\n", N, loop); Vnm_tstart(90, "VV Vset loop"); itotal = 0; for (j=0; j<loop; j++) for (vx=Gem_firstVV(thee);vx!=VNULL;vx=Gem_nextVV(thee)) itotal += ( j * VV_id(vx) * (j + VV_id(vx)) ); Vnm_tstop(90, "VV Vset loop"); Vnm_print(0,"Gem_speedChk: Freeing normal double vec...\n"); Vnm_tstart(90, "Free double vec"); Vmem_free( thee->vmem, N, sizeof(double), (void**)&d0 ); Vnm_tstop(90, "Free double vec"); Vnm_print(0,"Gem_speedChk: Freeing normal int vec...\n"); Vnm_tstart(90, "Free int vec"); Vmem_free( thee->vmem, N, sizeof(int), (void**)&i0 ); Vnm_tstop(90, "Free int vec"); Vnm_print(0,"Gem_speedChk: Freeing normal VV vec...\n"); Vnm_tstart(90, "Free VV vec"); Vmem_free( thee->vmem, N, sizeof(VV), (void**)&v0 ); Vnm_tstop(90, "Free VV vec"); Vnm_print(0,"Gem_speedChk: Freeing Vset VV vec...\n"); Vnm_tstart(90, "Free VV Vset"); for (i=0; i<N; i++) Gem_destroyVV(thee); Vnm_tstop(90, "Free VV Vset"); }