/* * *************************************************************************** * Routine: Bnode_ctor * * Purpose: The Block node constructor. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC Bnode* Bnode_ctor(Vmem *vmem, int pnumB, int pnumR[MAXV]) { int i; Bnode *thee = VNULL; VDEBUGIO("Bnode_ctor: CREATING object.."); thee = Vmem_malloc( VNULL, 1, sizeof(Bnode) ); if (vmem == VNULL) { thee->vmem = Vmem_ctor( "Bnode" ); thee->iMadeVmem = 1; } else { thee->vmem = vmem; thee->iMadeVmem = 0; } VDEBUGIO("..done.\n"); thee->numB = pnumB; for (i=0; i<thee->numB; i++) { thee->ND[i] = Node_ctor(vmem, pnumR[i]); } return thee; }
/* * *************************************************************************** * Routine: MCsh_ctor * * Purpose: The MCsh constructor. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC MCsh* MCsh_ctor(PDE *tpde, int argc, char **argv) { MCsh *thee = VNULL; VDEBUGIO("MCsh_ctor: CREATING object.."); thee = Vmem_malloc( VNULL, 1, sizeof(MCsh) ); thee->vmem = Vmem_ctor( "MCsh" ); VDEBUGIO("..done.\n"); /* the environment and shell object */ thee->PR[0] = '\0'; thee->USER_shell = VNULL; thee->vsh = Vsh_ctor(VNULL, argc, argv); /* the differential equation object */ thee->pde = tpde; /* the geometry manager object */ thee->gm = Gem_ctor( VNULL, thee->pde ); /* the approximation object */ thee->aprx = Aprx_ctor( VNULL, thee->gm, thee->pde ); /* the algebra manager object */ thee->am = AM_ctor( VNULL, thee->aprx ); return thee; }
/* * *************************************************************************** * Routine: Vmem_realloc * * Purpose: A logged version of realloc (using this is usually a bad idea). * * Author: Michael Holst * *************************************************************************** */ VPUBLIC void *Vmem_realloc(Vmem *thee, size_t num, size_t size, void **ram, size_t newNum) { void *tee = Vmem_malloc(thee, newNum, size); memcpy(tee, (*ram), size*VMIN2(num,newNum)); Vmem_free(thee, num, size, ram); return tee; }
/* /////////////////////////////////////////////////////////////////////////// // Routine: Vopot_ctor // Author: Nathan Baker /////////////////////////////////////////////////////////////////////////// */ VPUBLIC Vopot* Vopot_ctor(Vmgrid *mgrid, Vpbe *pbe, Vbcfl bcfl) { Vopot *thee = VNULL; thee = Vmem_malloc(VNULL, 1, sizeof(Vopot)); VASSERT(thee != VNULL); VASSERT(Vopot_ctor2(thee, mgrid, pbe, bcfl)); return thee; }
VPUBLIC Vcsm* Vcsm_ctor(Valist *alist, Gem *gm) { /* Set up the structure */ Vcsm *thee = VNULL; thee = Vmem_malloc(VNULL, 1, sizeof(Vcsm) ); VASSERT( thee != VNULL); VASSERT( Vcsm_ctor2(thee, alist, gm)); return thee; }
VPUBLIC FEMparm* FEMparm_ctor(FEMparm_CalcType type) { /* Set up the structure */ FEMparm *thee = VNULL; thee = Vmem_malloc(VNULL, 1, sizeof(FEMparm)); VASSERT( thee != VNULL); VASSERT( FEMparm_ctor2(thee, type) ); return thee; }
VPUBLIC BEMparm* BEMparm_ctor(BEMparm_CalcType type) { /* Set up the structure */ BEMparm *thee = VNULL; thee = (BEMparm*)Vmem_malloc(VNULL, 1, sizeof(BEMparm)); VASSERT( thee != VNULL); VASSERT( BEMparm_ctor2(thee, type) == VRC_SUCCESS ); return thee; }
VPUBLIC Vgreen* Vgreen_ctor(Valist *alist) { /* Set up the structure */ Vgreen *thee = VNULL; thee = (Vgreen *)Vmem_malloc(VNULL, 1, sizeof(Vgreen) ); VASSERT( thee != VNULL); VASSERT( Vgreen_ctor2(thee, alist)); return thee; }
VPUBLIC VclistCell* VclistCell_ctor(int natoms) { VclistCell *thee = VNULL; /* Set up the structure */ thee = Vmem_malloc(VNULL, 1, sizeof(VclistCell)); VASSERT( thee != VNULL); VASSERT( VclistCell_ctor2(thee, natoms) == VRC_SUCCESS ); return thee; }
/* * *************************************************************************** * Routine: Vcom_ctor * * Purpose: Construct the communications object * * Notes: This routine sets up data members of class and initializes MPI. * * Author: Nathan Baker and Michael Holst * *************************************************************************** */ VPUBLIC Vcom* Vcom_ctor(int commtype) { int rc; Vcom *thee = VNULL; /* Set up the structure */ thee = Vmem_malloc( VNULL, 1, sizeof(Vcom) ); thee->core = Vmem_malloc( VNULL, 1, sizeof(Vcom_core) ); /* Call the real constructor */ rc = Vcom_ctor2(thee, commtype); /* Destroy the guy if something went wrong */ if (rc == 0) { Vmem_free( VNULL, 1, sizeof(Vcom_core), (void**)&(thee->core) ); Vmem_free( VNULL, 1, sizeof(Vcom), (void**)&thee ); } return thee; }
/* Main (FORTRAN stub) constructor */ VPUBLIC Vrc_Codes Vclist_ctor2(Vclist *thee, Valist *alist, double max_radius, int npts[VAPBS_DIM], Vclist_DomainMode mode, double lower_corner[VAPBS_DIM], double upper_corner[VAPBS_DIM]) { int i; VclistCell *cell; /* Check and store parameters */ if ( Vclist_storeParms(thee, alist, max_radius, npts, mode, lower_corner, upper_corner) == VRC_FAILURE ) { Vnm_print(2, "Vclist_ctor2: parameter check failed!\n"); return VRC_FAILURE; } /* Set up memory */ thee->vmem = Vmem_ctor("APBS::VCLIST"); if (thee->vmem == VNULL) { Vnm_print(2, "Vclist_ctor2: memory object setup failed!\n"); return VRC_FAILURE; } /* Set up cells */ thee->cells = Vmem_malloc( thee->vmem, thee->n, sizeof(VclistCell) ); if (thee->cells == VNULL) { Vnm_print(2, "Vclist_ctor2: Failed allocating %d VclistCell objects!\n", thee->n); return VRC_FAILURE; } for (i=0; i<thee->n; i++) { cell = &(thee->cells[i]); cell->natoms = 0; } /* Set up the grid */ if ( Vclist_setupGrid(thee) == VRC_FAILURE ) { Vnm_print(2, "Vclist_ctor2: grid setup failed!\n"); return VRC_FAILURE; } /* Assign atoms to grid cells */ if (Vclist_assignAtoms(thee) == VRC_FAILURE) { Vnm_print(2, "Vclist_ctor2: atom assignment failed!\n"); return VRC_FAILURE; } return VRC_SUCCESS; }
VPUBLIC Vclist* Vclist_ctor(Valist *alist, double max_radius, int npts[VAPBS_DIM], Vclist_DomainMode mode, double lower_corner[VAPBS_DIM], double upper_corner[VAPBS_DIM]) { Vclist *thee = VNULL; /* Set up the structure */ thee = Vmem_malloc(VNULL, 1, sizeof(Vclist) ); VASSERT( thee != VNULL); VASSERT( Vclist_ctor2(thee, alist, max_radius, npts, mode, lower_corner, upper_corner) == VRC_SUCCESS ); return thee; }
VPUBLIC int Vgreen_coulomb(Vgreen *thee, int npos, double *x, double *y, double *z, double *val) { Vatom *atom; double *apos, charge, dist, dx, dy, dz, scale; double *q, qtemp, fx, fy, fz; int iatom, ipos; if (thee == VNULL) { Vnm_print(2, "Vgreen_coulomb: Got NULL thee!\n"); return 0; } for (ipos=0; ipos<npos; ipos++) val[ipos] = 0.0; #ifdef HAVE_TREE /* Allocate charge array (if necessary) */ if (Valist_getNumberAtoms(thee->alist) > 1) { if (npos > 1) { q = VNULL; q = Vmem_malloc(thee->vmem, npos, sizeof(double)); if (q == VNULL) { Vnm_print(2, "Vgreen_coulomb: Error allocating charge array!\n"); return 0; } } else { q = &(qtemp); } for (ipos=0; ipos<npos; ipos++) q[ipos] = 1.0; /* Calculate */ treecalc(thee, x, y, z, q, npos, val, thee->xp, thee->yp, thee->zp, thee->qp, thee->np, &fx, &fy, &fz, 1, 1, thee->np); } else return Vgreen_coulomb_direct(thee, npos, x, y, z, val); /* De-allocate charge array (if necessary) */ if (npos > 1) Vmem_free(thee->vmem, npos, sizeof(double), (void **)&q); scale = Vunit_ec/(4*Vunit_pi*Vunit_eps0*1.0e-10); for (ipos=0; ipos<npos; ipos++) val[ipos] = val[ipos]*scale; return 1; #else /* ifdef HAVE_TREE */ return Vgreen_coulomb_direct(thee, npos, x, y, z, val); #endif }
/* * *************************************************************************** * Routine: Bmat_ctor * * Purpose: The block sparse matrix constructor. * * Notes: This constructor only does minimal initialization of a Bmat * object after creating the object itself; it doesn't create * any storage for either the integer structure arrays or the * nonzero storage arrays. * * This constructor only fixes the number of blocks and the * numbers of rows and columns in each block; the nonzero * structure is not set yet. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC Bmat* Bmat_ctor(Vmem *vmem, const char *name, int pnumB, int pnumR[MAXV], int pnumC[MAXV], MATmirror pmirror[MAXV][MAXV]) { int p,q; char bname[15]; Bmat *thee = VNULL; VDEBUGIO("Bmat_ctor: CREATING object.."); thee = Vmem_malloc( VNULL, 1, sizeof(Bmat) ); if (vmem == VNULL) { thee->vmem = Vmem_ctor( "Bmat" ); thee->iMadeVmem = 1; } else { thee->vmem = vmem; thee->iMadeVmem = 0; } strncpy(thee->name, name, 10); thee->coarse = VNULL; thee->fine = VNULL; thee->numB = pnumB; for (p=0; p<thee->numB; p++) { for (q=0; q<thee->numB; q++) { thee->mirror[p][q] = pmirror[p][q]; thee->AD[p][q] = VNULL; } } /* create upper-triangular blocks before the lower-triangular blocks */ for (p=0; p<thee->numB; p++) { for (q=0; q<thee->numB; q++) { if ( !thee->mirror[p][q] ) { sprintf(bname, "%s_%d%d", thee->name, p, q); thee->AD[p][q] = Mat_ctor(thee->vmem,bname,pnumR[p],pnumC[q]); } else { /* ( thee->mirror[p][q] ) */ /* first make sure that the mirror of this guy exists! */ VASSERT( !thee->mirror[q][p] ); VASSERT( thee->AD[q][p] != VNULL ); thee->AD[p][q] = thee->AD[q][p]; } } } VDEBUGIO("..done.\n"); return thee; }
/* * *************************************************************************** * Routine: Vmem_ctor * * Purpose: Construct the dynamic memory allocation logging object. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC Vmem *Vmem_ctor(char *name) { Vmem *thee; thee = Vmem_malloc( VNULL, 1, sizeof(Vmem) ); VASSERT( thee != VNULL ); strncpy( thee->name, name, VMAX_ARGLEN ); thee->mallocBytes = 0; thee->freeBytes = 0; thee->highWater = 0; thee->mallocAreas = 0; return thee; }
/* /////////////////////////////////////////////////////////////////////////// // Routine: Vpee_ctor // // Author: Nathan Baker /////////////////////////////////////////////////////////////////////////// */ VPUBLIC Vpee* Vpee_ctor(Gem *gm, int localPartID, int killFlag, double killParam ) { Vpee *thee = VNULL; /* Set up the structure */ thee = Vmem_malloc(VNULL, 1, sizeof(Vpee) ); VASSERT( thee != VNULL); VASSERT( Vpee_ctor2(thee, gm, localPartID, killFlag, killParam)); return thee; }
/* /////////////////////////////////////////////////////////////////////////// // Routine: PDE_ctor // // Purpose: Construct the differential equation object. // // Speed: This function is called by MC one time during setup, // and does not have to be particularly fast. // // Author: Michael Holst /////////////////////////////////////////////////////////////////////////// */ VPUBLIC PDE* myPDE_ctor(void) { int i; PDE *thee = VNULL; /* create some space for the pde object */ thee = Vmem_malloc( VNULL, 1, sizeof(PDE) ); /* PDE-specific parameters and function pointers */ thee->initAssemble = initAssemble; /* once-per-assembly initialization */ thee->initElement = initElement; /* once-per-element initialization */ thee->initFace = initFace; /* once-per-face initialization */ thee->initPoint = initPoint; /* once-per-point initialization */ thee->Fu = Fu; /* nonlinear strong form */ thee->Ju = Ju; /* nonlinear energy functional */ thee->Fu_v = Fu_v; /* nonlinear weak form */ thee->DFu_wv = DFu_wv; /* bilinear linearization weak form */ thee->p_wv = p_wv; /* bilinear mass density form */ thee->delta = delta; /* delta function source term */ thee->u_D = u_D; /* dirichlet func and initial guess */ thee->u_T = u_T; /* analytical soln for testing */ thee->vec = 1; /* unknowns per spatial point; */ thee->sym[0][0] = 0; /* symmetries of bilinear form */ thee->est[0] = 1.0; /* error estimator weights */ for (i=0; i<VMAX_BDTYPE; i++) /* boundary type remappings */ thee->bmap[0][i] = i; /* Manifold-specific function pointers */ thee->bisectEdge = bisectEdge; /* edge bisection rule */ thee->mapBoundary = mapBoundary; /* boundary recovery rule */ thee->markSimplex = markSimplex; /* simplex marking rule */ thee->oneChart = oneChart; /* coordinate transformations */ /* Element-specific function pointers */ thee->simplexBasisInit = simplexBasisInit; /* initialization of bases */ thee->simplexBasisForm = simplexBasisForm; /* form trial & test bases */ /* Special DYN structure */ PDE_initDyn( thee ); /* return the new pde object */ return thee; }
/* * *************************************************************************** * Routine: Vmpi_ctor * * Purpose: The Vmpi constructor. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC Vmpi* Vmpi_ctor(void) { #if defined(HAVE_MPI_H) int dummy; #endif Vmpi *thee = VNULL; VDEBUGIO("Vmpi_ctor: CREATING object.."); thee = Vmem_malloc( VNULL, 1, sizeof(Vmpi) ); thee->mpi_rank = 0; thee->mpi_size = 0; #if defined(HAVE_MPI_H) VASSERT( MPI_SUCCESS == MPI_Initialized(&dummy) ); VASSERT( MPI_SUCCESS == MPI_Comm_rank(MPI_COMM_WORLD, &(thee->mpi_rank)) ); VASSERT( MPI_SUCCESS == MPI_Comm_size(MPI_COMM_WORLD, &(thee->mpi_size)) ); Vnm_setIoTag(thee->mpi_rank, thee->mpi_size); Vnm_print(2,"Vmpi_ctor: process %d of %d is ALIVE!\n", thee->mpi_rank, thee->mpi_size); #endif VDEBUGIO("..done.\n"); return thee; }
VPUBLIC Vrc_Codes VclistCell_ctor2(VclistCell *thee, int natoms) { if (thee == VNULL) { Vnm_print(2, "VclistCell_ctor2: NULL thee!\n"); return VRC_FAILURE; } thee->natoms = natoms; if (thee->natoms > 0) { thee->atoms = Vmem_malloc(VNULL, natoms, sizeof(Vatom *)); if (thee->atoms == VNULL) { Vnm_print(2, "VclistCell_ctor2: unable to allocate space for %d atom pointers!\n", natoms); return VRC_FAILURE; } } return VRC_SUCCESS; }
/* * *************************************************************************** * Routine: Gem_ctor * * Purpose: Geometry manager constructor. * * Author: Michael Holst * *************************************************************************** */ VPUBLIC Gem* Gem_ctor(Vmem *vmem, PDE *tpde) { int i; Gem* thee = VNULL; thee = Vmem_malloc( VNULL, 1, sizeof(Gem) ); if (vmem == VNULL) { thee->vmem = Vmem_ctor( "Gem" ); thee->iMadeVmem = 1; } else { thee->vmem = vmem; thee->iMadeVmem = 0; } VDEBUGIO("Gem_ctor: CREATING object.."); VDEBUGIO("..done.\n"); thee->xUp = defaultxUpFunc; thee->xUpFlag = 0; if (tpde != VNULL) { thee->iMadePDE = 0; thee->pde = tpde; } else { thee->iMadePDE = 1; thee->pde = PDE_ctor_default(); } thee->vertices = VNULL; thee->edges = VNULL; thee->simplices = VNULL; for (i=0; i<VMAXSQ; i++) { thee->sQueM[i] = VNULL; } Gem_reset(thee, 0, 0); return thee; }
/* * *************************************************************************** * Routine: Slu_ctor * * Purpose: The Slu constructor. * * Author: Michael Holst and Stephen Bond * *************************************************************************** */ VPUBLIC Slu* Slu_ctor(Vmem *vmem, int skey, int m, int n, int nnz, int *ia, int *ja, double *a) { Slu* thee = VNULL; VDEBUGIO("Slu_ctor: CREATING object.."); thee = Vmem_malloc( VNULL, 1, sizeof(Slu) ); if (vmem == VNULL) { thee->vmem = Vmem_ctor( "Slu" ); thee->iMadeVmem = 1; } else { thee->vmem = vmem; thee->iMadeVmem = 0; } /* dimensions */ thee->m = m; thee->n = n; /* storage key */ thee->skey = skey; /* pointer to factors */ thee->work = VNULL; /* arrays */ thee->a = a; thee->ia = ia; thee->ja = ja; VDEBUGIO("..done.\n"); /* return */ return thee; }
void apbsnlinduce_(double uinp[maxatm][3], double fld[maxatm][3]){ /* Misc. pointers to APBS data structures */ Vpmg *pmg[NOSH_MAXCALC]; Vpmgp *pmgp[NOSH_MAXCALC]; Vpbe *pbe[NOSH_MAXCALC]; MGparm *mgparm = VNULL; PBEparm *pbeparm = VNULL; Vatom *atom = VNULL; /* Observables and unit conversion */ double field[3]; double sign,kT,electric; /* Potential Vgrid construction */ double nx,ny,nz,hx,hy,hzed,xmin,ymin,zmin; double *data; /* Loop variables */ int i,j; VASSERT(nosh != VNULL); for (i=0; i<NOSH_MAXCALC; i++) { pmg[i] = VNULL; pmgp[i] = VNULL; pbe[i] = VNULL; } /* Read TINKER induce input data into Vatom instances. */ for (i=0; i < alist[0]->number; i++){ atom = Valist_getAtom(alist[0],i); Vatom_setNLInducedDipole(atom, uinp[i]); for (j=0;j<3;j++){ fld[i][j] = 0.0; } } /* Solve the LPBE for the homogeneous system, then solvated. */ for (i=0; i<2; i++) { pmg[i] = VNULL; pmgp[i] = VNULL; pbe[i] = VNULL; /* Useful local variables */ mgparm = nosh->calc[i]->mgparm; pbeparm = nosh->calc[i]->pbeparm; if (!MGparm_check(mgparm)){ printf("MGparm Check failed\n"); exit(-1); } if (!PBEparm_check(pbeparm)){ printf("PBEparm Check failed\n"); exit(-1); } /* Set up problem */ mgparm->chgs = VCM_NLINDUCED; if (!initMG(i, nosh, mgparm, pbeparm, realCenter, pbe, alist, dielXMap, dielYMap, dielZMap, kappaMap, chargeMap, pmgp, pmg, potMap)) { Vnm_tprint( 2, "Error setting up MG calculation!\n"); return; } /* Solve the PDE */ if (solveMG(nosh, pmg[i], mgparm->type) != 1) { Vnm_tprint(2, "Error solving PDE!\n"); return; } /* Set partition information for observables and I/O */ if (setPartMG(nosh, mgparm, pmg[i]) != 1) { Vnm_tprint(2, "Error setting partition info!\n"); return; } /* Save the potential due to non-local induced dipoles */ nx = pmg[i]->pmgp->nx; ny = pmg[i]->pmgp->ny; nz = pmg[i]->pmgp->nz; hx = pmg[i]->pmgp->hx; hy = pmg[i]->pmgp->hy; hzed = pmg[i]->pmgp->hzed; xmin = pmg[i]->pmgp->xmin; ymin = pmg[i]->pmgp->ymin; zmin = pmg[i]->pmgp->zmin; if (nlIndU[i] == VNULL) { data = Vmem_malloc(VNULL, nx*ny*nz, sizeof(double)); Vpmg_fillArray(pmg[i], data, VDT_POT, 0.0, pbeparm->pbetype, pbeparm); nlIndU[i] = Vgrid_ctor(nx,ny,nz,hx,hy,hzed,xmin,ymin,zmin,data); nlIndU[i]->readdata = 1; // set readata flag to have dtor free data } else { data = nlIndU[i]->data; Vpmg_fillArray(pmg[i], data, VDT_POT, 0.0, pbeparm->pbetype, pbeparm); } if (i == 0){ sign = -1.0; } else { sign = 1.0; } for (j=0; j < alist[0]->number; j++){ Vpmg_fieldSpline4(pmg[i], j, field); fld[j][0] += sign * field[0]; fld[j][1] += sign * field[1]; fld[j][2] += sign * field[2]; } } /* load results into the return arrays in electron**2/Angstrom /* kT in kcal/mol */ kT = Vunit_kb * (1e-3) * Vunit_Na * 298.15 / 4.184; // electric: conversion from electron**2/Angstrom to Kcal/mol electric = 332.063709; for (i=0; i<alist[0]->number; i++){ fld[i][0] *= kT / electric; fld[i][1] *= kT / electric; fld[i][2] *= kT / electric; } killMG(nosh, pbe, pmgp, pmg); }
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; }
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: 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: 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; }
VPRIVATE int treesetup(Vgreen *thee) { #ifdef HAVE_TREE double dist_tol = FMM_DIST_TOL; int iflag = FMM_IFLAG; double order = FMM_ORDER; int theta = FMM_THETA; int shrink = FMM_SHRINK; int maxparnode = FMM_MAXPARNODE; int minlevel = FMM_MINLEVEL; int maxlevel = FMM_MAXLEVEL; int level = 0; int one = 1; Vatom *atom; double xyzminmax[6], *pos; int i; /* Set up particle arrays with atomic coordinates and charges */ Vnm_print(0, "treesetup: Initializing FMM particle arrays...\n"); thee->np = Valist_getNumberAtoms(thee->alist); thee->xp = VNULL; thee->xp = (double *)Vmem_malloc(thee->vmem, thee->np, sizeof(double)); if (thee->xp == VNULL) { Vnm_print(2, "Vgreen_ctor2: Failed to allocate %d*sizeof(double)!\n", thee->np); return 0; } thee->yp = VNULL; thee->yp = (double *)Vmem_malloc(thee->vmem, thee->np, sizeof(double)); if (thee->yp == VNULL) { Vnm_print(2, "Vgreen_ctor2: Failed to allocate %d*sizeof(double)!\n", thee->np); return 0; } thee->zp = VNULL; thee->zp = (double *)Vmem_malloc(thee->vmem, thee->np, sizeof(double)); if (thee->zp == VNULL) { Vnm_print(2, "Vgreen_ctor2: Failed to allocate %d*sizeof(double)!\n", thee->np); return 0; } thee->qp = VNULL; thee->qp = (double *)Vmem_malloc(thee->vmem, thee->np, sizeof(double)); if (thee->qp == VNULL) { Vnm_print(2, "Vgreen_ctor2: Failed to allocate %d*sizeof(double)!\n", thee->np); return 0; } for (i=0; i<thee->np; i++) { atom = Valist_getAtom(thee->alist, i); pos = Vatom_getPosition(atom); thee->xp[i] = pos[0]; thee->yp[i] = pos[1]; thee->zp[i] = pos[2]; thee->qp[i] = Vatom_getCharge(atom); } Vnm_print(0, "treesetup: Setting things up...\n"); F77SETUP(thee->xp, thee->yp, thee->zp, &(thee->np), &order, &theta, &iflag, &dist_tol, xyzminmax, &(thee->np)); Vnm_print(0, "treesetup: Initializing levels...\n"); F77INITLEVELS(&minlevel, &maxlevel); Vnm_print(0, "treesetup: Creating tree...\n"); F77CREATE_TREE(&one, &(thee->np), thee->xp, thee->yp, thee->zp, thee->qp, &shrink, &maxparnode, xyzminmax, &level, &(thee->np)); return 1; #else /* ifdef HAVE_TREE */ Vnm_print(2, "treesetup: Error! APBS not linked with treecode!\n"); return 0; #endif /* ifdef HAVE_TREE */ }
VPUBLIC int Vgreen_coulombD(Vgreen *thee, int npos, double *x, double *y, double *z, double *pot, double *gradx, double *grady, double *gradz) { Vatom *atom; double *apos, charge, dist, dist2, idist3, dy, dz, dx, scale; double *q, qtemp; int iatom, ipos; if (thee == VNULL) { Vnm_print(2, "Vgreen_coulombD: Got VNULL thee!\n"); return 0; } for (ipos=0; ipos<npos; ipos++) { pot[ipos] = 0.0; gradx[ipos] = 0.0; grady[ipos] = 0.0; gradz[ipos] = 0.0; } #ifdef HAVE_TREE if (Valist_getNumberAtoms(thee->alist) > 1) { if (npos > 1) { q = VNULL; q = Vmem_malloc(thee->vmem, npos, sizeof(double)); if (q == VNULL) { Vnm_print(2, "Vgreen_coulomb: Error allocating charge array!\n"); return 0; } } else { q = &(qtemp); } for (ipos=0; ipos<npos; ipos++) q[ipos] = 1.0; /* Calculate */ treecalc(thee, x, y, z, q, npos, pot, thee->xp, thee->yp, thee->zp, thee->qp, thee->np, gradx, grady, gradz, 2, npos, thee->np); /* De-allocate charge array (if necessary) */ if (npos > 1) Vmem_free(thee->vmem, npos, sizeof(double), (void **)&q); } else return Vgreen_coulombD_direct(thee, npos, x, y, z, pot, gradx, grady, gradz); scale = Vunit_ec/(4*VPI*Vunit_eps0*(1.0e-10)); for (ipos=0; ipos<npos; ipos++) { gradx[ipos] = gradx[ipos]*scale; grady[ipos] = grady[ipos]*scale; gradz[ipos] = gradz[ipos]*scale; pot[ipos] = pot[ipos]*scale; } return 1; #else /* ifdef HAVE_TREE */ return Vgreen_coulombD_direct(thee, npos, x, y, z, pot, gradx, grady, gradz); #endif }
/* * *************************************************************************** * Routine: Mat_initStructureLN * * Purpose: Initialize the nonzero structure given structure information. * * Notes: This routine actually does the storage creation for all * internal Vset, link arrays, and link pointer arrays. * * Author: Stephen Bond and Michael Holst * *************************************************************************** */ VPUBLIC void Mat_initStructureLN(Mat *thee, MATformat frmt, MATsym sym) { int i, sizeA; LinkA *mt; LinkRCS *mtX; /* use an upperbound for the maximum size of the Vset */ sizeA = (thee->numR)*LN_MAX_ENTRIES_PER_ROW; /* build framework from input */ thee->format = frmt; thee->sym = sym; /* initialize the state */ thee->state = ZERO_STATE; /* Format-dependent consistency checks and setup */ switch (thee->format) { /* RLN consistency checks and setup */ case RLN_FORMAT: /* RLN valid for nonsymmetry */ VASSERT( thee->sym==ISNOT_SYM ); thee->lnkU = Vset_ctor( thee->vmem, "U", sizeof(LinkA), sizeA, 0 ); /* create an empty entry to start each row of the matrix */ for (i=0; i<(thee->numR); i++) { mt=(LinkA*)Vset_create(thee->lnkU); mt->idx = -1; mt->val = 0.; mt->next = VNULL; } thee->numO = 0; break; /* CLN consistency checks and setup */ case CLN_FORMAT: /* RLN valid for nonsymmetry */ VASSERT( thee->sym==ISNOT_SYM ); thee->lnkL = Vset_ctor( thee->vmem, "L", sizeof(LinkA), sizeA, 0 ); /* create an empty entry to start each column of the matrix */ for (i=0; i<(thee->numC); i++) { mt=(LinkA*)Vset_create(thee->lnkL); mt->idx = -1; mt->val = 0.; mt->next = VNULL; } thee->numO = 0; break; /* XLN consistency checks and setup */ case XLN_FORMAT: switch( sym ) { case IS_SYM: thee->lnkU = Vset_ctor( thee->vmem,"X",sizeof(LinkRCS),sizeA,0 ); thee->numA = thee->numR; thee->xln = Vmem_malloc( thee->vmem,thee->numA,sizeof(LinkRCS) ); /* create a zero entries to start the diagonal of the matrix */ for (i=0; i<(thee->numR); i++) { mtX = &( ((LinkRCS*) (thee->xln))[i] ); mtX->idxT = i; mtX->idx = i; mtX->nxtT = VNULL; mtX->next = VNULL; mtX->val = 0.; } break; case STRUC_SYM: thee->lnkU = Vset_ctor( thee->vmem,"X",sizeof(LinkRCN),sizeA,0 ); thee->numA = thee->numR; thee->xln = Vmem_malloc( thee->vmem,thee->numA,sizeof(LinkRCS) ); /* create a zero entries to start the diagonal of the matrix */ for (i=0; i<(thee->numR); i++) { mtX = &( ((LinkRCS*) (thee->xln))[i] ); mtX->idxT = i; mtX->idx = i; mtX->nxtT = VNULL; mtX->next = VNULL; mtX->val = 0.; } break; case ISNOT_SYM: thee->lnkU = Vset_ctor( thee->vmem,"X",sizeof(LinkRCS),sizeA,0 ); thee->numA = thee->numR + thee->numC; thee->xln = Vmem_malloc( thee->vmem,thee->numA,sizeof(LinkRC*) ); thee->xlnt = (void*) &( ((LinkRC**) (thee->xln))[thee->numR] ); break; default: VASSERT(0); break; } thee->numO = 0; break; default: VASSERT(0); break; } }
void apbsempole_(int *natom, double x[maxatm][3], double rad[maxatm], double rpole[maxatm][13], double *total, double energy[maxatm], double fld[maxatm][3], double rff[maxatm][3], double rft[maxatm][3]) { /* Misc. pointers to APBS data structures */ Vpmg *pmg[NOSH_MAXCALC]; Vpmgp *pmgp[NOSH_MAXCALC]; Vpbe *pbe[NOSH_MAXCALC]; MGparm *mgparm = VNULL; PBEparm *pbeparm = VNULL; Vatom *atom = VNULL; /* Vgrid configuration for the kappa and dielectric maps */ double nx,ny,nz,hx,hy,hzed,xmin,ymin,zmin; double *data; double zkappa2, epsp, epsw; /* Loop indeces */ int i,j; /* Observables and unit conversion */ double sign, force[3], torque[3], field[3]; double kT,electric,debye; double charge, dipole[3], quad[9]; debye = 4.8033324; for (i=0; i<NOSH_MAXCALC; i++) { pmg[i] = VNULL; pmgp[i] = VNULL; pbe[i] = VNULL; } /* Kill the saved potential Vgrids */ for (i=0; i<2; i++){ if (permU[i] != VNULL) Vgrid_dtor(&permU[i]); if (indU[i] != VNULL) Vgrid_dtor(&indU[i]); if (nlIndU[i] != VNULL) Vgrid_dtor(&nlIndU[i]); } /* Kill the old atom list */ if (alist[0] != VNULL) { Valist_dtor(&alist[0]); } /* Create a new atom list (mol == 1) */ if (alist[0] == VNULL) { alist[0] = Valist_ctor(); alist[0]->atoms = Vmem_malloc(alist[0]->vmem, *natom, (sizeof(Vatom))); alist[0]->number = *natom; } /* Read TINKER input data into Vatom instances. */ for (i=0; i < alist[0]->number; i++){ atom = Valist_getAtom(alist[0],i); Vatom_setAtomID(atom, i); Vatom_setPosition(atom, x[i]); Vatom_setRadius(atom, rad[i]); charge = rpole[i][0]; Vatom_setCharge(atom, charge); dipole[0] = rpole[i][1]; dipole[1] = rpole[i][2]; dipole[2] = rpole[i][3]; Vatom_setDipole(atom, dipole); quad[0] = rpole[i][4]; quad[1] = rpole[i][5]; quad[2] = rpole[i][6]; quad[3] = rpole[i][7]; quad[4] = rpole[i][8]; quad[5] = rpole[i][9]; quad[6] = rpole[i][10]; quad[7] = rpole[i][11]; quad[8] = rpole[i][12]; Vatom_setQuadrupole(atom, quad); /* Useful check printf(" %i %f (%f,%f,%f)\n",i,rad[i], x[i][0], x[i][1], x[i][2]); printf(" %f\n %f,%f,%f\n", charge, dipole[0], dipole[1], dipole[2]); printf(" %f\n", quad[0]); printf(" %f %f\n", quad[3], quad[4]); printf(" %f %f %f\n", quad[6], quad[7], quad[8]); */ energy[i] = 0.0; for (j=0;j<3;j++){ fld[i][j] = 0.0; rff[i][j] = 0.0; rft[i][j] = 0.0; } } nosh->nmol = 1; Valist_getStatistics(alist[0]); /* Only call the setupCalc routine once, so that we can reuse this nosh object */ if (nosh->ncalc < 2) { if (NOsh_setupElecCalc(nosh, alist) != 1) { printf("Error setting up calculations\n"); exit(-1); } } /* Solve the LPBE for the homogeneous and then solvated states */ for (i=0; i<2; i++) { /* Useful local variables */ mgparm = nosh->calc[i]->mgparm; pbeparm = nosh->calc[i]->pbeparm; /* Just to be robust */ if (!MGparm_check(mgparm)){ printf("MGparm Check failed\n"); printMGPARM(mgparm, realCenter); exit(-1); } if (!PBEparm_check(pbeparm)){ printf("PBEparm Check failed\n"); printPBEPARM(pbeparm); exit(-1); } /* Set up the problem */ mgparm->chgs = VCM_PERMANENT; if (!initMG(i, nosh, mgparm, pbeparm, realCenter, pbe, alist, dielXMap, dielYMap, dielZMap, kappaMap, chargeMap, pmgp, pmg, potMap)) { Vnm_tprint( 2, "Error setting up MG calculation!\n"); return; } /* Solve the PDE */ if (solveMG(nosh, pmg[i], mgparm->type) != 1) { Vnm_tprint(2, "Error solving PDE!\n"); return; } /* Set partition information for observables and I/O */ /* Note - parallel operation has NOT been tested. */ if (setPartMG(nosh, mgparm, pmg[i]) != 1) { Vnm_tprint(2, "Error setting partition info!\n"); return; } nx = pmg[i]->pmgp->nx; ny = pmg[i]->pmgp->ny; nz = pmg[i]->pmgp->nz; hx = pmg[i]->pmgp->hx; hy = pmg[i]->pmgp->hy; hzed = pmg[i]->pmgp->hzed; xmin = pmg[i]->pmgp->xmin; ymin = pmg[i]->pmgp->ymin; zmin = pmg[i]->pmgp->zmin; /* Save dielectric/kappa maps into Vgrids, then change the nosh * data structure to think it read these maps in from a file. * The goal is to save setup time during convergence of the * induced dipoles. This is under consideration... * */ /* // X (shifted) data = Vmem_malloc(mem, nx*ny*nz, sizeof(double)); Vpmg_fillArray(pmg[i], data, VDT_DIELX, 0.0, pbeparm->pbetype); dielXMap[i] = Vgrid_ctor(nx,ny,nz,hx,hy,hzed, xmin + 0.5*hx,ymin,zmin,data); dielXMap[i]->readdata = 1; // Y (shifted) data = Vmem_malloc(mem, nx*ny*nz, sizeof(double)); Vpmg_fillArray(pmg[i], data, VDT_DIELY, 0.0, pbeparm->pbetype); dielYMap[i] = Vgrid_ctor(nx,ny,nz,hx,hy,hzed, xmin,ymin + 0.5*hy,zmin,data); dielYMap[i]->readdata = 1; // Z (shifted) data = Vmem_malloc(mem, nx*ny*nz, sizeof(double)); Vpmg_fillArray(pmg[i], data, VDT_DIELZ, 0.0, pbeparm->pbetype); dielZMap[i] = Vgrid_ctor(nx,ny,nz,hx,hy,hzed, xmin,ymin,zmin + 0.5*hzed,data); dielZMap[i]->readdata = 1; // Kappa data = Vmem_malloc(mem, nx*ny*nz, sizeof(double)); Vpmg_fillArray(pmg[i], data, VDT_KAPPA, 0.0, pbeparm->pbetype); kappaMap[i] = Vgrid_ctor(nx,ny,nz,hx,hy,hzed,xmin,ymin,zmin,data); kappaMap[i]->readdata = 1; // Update the pbeparam structure, since we now have // dielectric and kappap maps pbeparm->useDielMap = 1; pbeparm->dielMapID = i + 1; pbeparm->useKappaMap = 1; pbeparm->kappaMapID = i + 1; */ data = Vmem_malloc(mem, nx*ny*nz, sizeof(double)); Vpmg_fillArray(pmg[i], data, VDT_POT, 0.0, pbeparm->pbetype, pbeparm); permU[i] = Vgrid_ctor(nx,ny,nz,hx,hy,hzed,xmin,ymin,zmin,data); permU[i]->readdata = 1; // set readdata flag to have the dtor to free data if (i == 0){ sign = -1.0; } else { sign = 1.0; } /* Calculate observables */ for (j=0; j < alist[0]->number; j++){ energy[j] += sign * Vpmg_qfPermanentMultipoleEnergy(pmg[i], j); Vpmg_fieldSpline4(pmg[i], j, field); fld[j][0] += sign * field[0]; fld[j][1] += sign * field[1]; fld[j][2] += sign * field[2]; } if (!pmg[i]->pmgp->nonlin && (pmg[i]->surfMeth == VSM_SPLINE || pmg[i]->surfMeth == VSM_SPLINE3 || pmg[i]->surfMeth == VSM_SPLINE4)) { for (j=0; j < alist[0]->number; j++){ Vpmg_qfPermanentMultipoleForce(pmg[i], j, force, torque); rff[j][0] += sign * force[0]; rff[j][1] += sign * force[1]; rff[j][2] += sign * force[2]; rft[j][0] += sign * torque[0]; rft[j][1] += sign * torque[1]; rft[j][2] += sign * torque[2]; } kT = Vunit_kb * (1e-3) * Vunit_Na * 298.15 * 1.0/4.184; epsp = Vpbe_getSoluteDiel(pmg[i]->pbe); epsw = Vpbe_getSolventDiel(pmg[i]->pbe); if (VABS(epsp-epsw) > VPMGSMALL) { for (j=0; j < alist[0]->number; j++){ Vpmg_dbPermanentMultipoleForce(pmg[i], j, force); rff[j][0] += sign * force[0]; rff[j][1] += sign * force[1]; rff[j][2] += sign * force[2]; } } zkappa2 = Vpbe_getZkappa2(pmg[i]->pbe); if (zkappa2 > VPMGSMALL) { for (j=0; j < alist[0]->number; j++) { Vpmg_ibPermanentMultipoleForce(pmg[i], j, force); rff[j][0] += sign * force[0]; rff[j][1] += sign * force[1]; rff[j][2] += sign * force[2]; } } } } //nosh->ndiel = 2; //nosh->nkappa = 2; /* printf("Energy (multipole) %f Kcal/mol\n", *energy); printf("Energy (volume) %f Kcal/mol\n", evol * 0.5 * kT); */ // Convert results into kcal/mol units kT = Vunit_kb * (1e-3) * Vunit_Na * 298.15 * 1.0/4.184; // Electric converts from electron**2/Angstrom to kcal/mol electric = 332.063709; *total = 0.0; for (i=0; i<alist[0]->number; i++){ /* starting with the field in KT/e/Ang^2 multiply by kcal/mol/KT the field is then divided by "electric" to convert to e/Ang^2 */ energy[i] *= 0.5 * kT; *total += energy[i]; fld[i][0] *= kT / electric; fld[i][1] *= kT / electric; fld[i][2] *= kT / electric; rff[i][0] *= kT; rff[i][1] *= kT; rff[i][2] *= kT; rft[i][0] *= kT; rft[i][1] *= kT; rft[i][2] *= kT; } killMG(nosh, pbe, pmgp, pmg); }