static real *mk_nbfp(t_idef *idef,bool bBHAM) { real *nbfp; int i,j,k,atnr; atnr=idef->atnr; if (bBHAM) { snew(nbfp,3*atnr*atnr); for(i=k=0; (i<atnr); i++) { for(j=0; (j<atnr); j++,k++) { BHAMA(nbfp,atnr,i,j) = idef->iparams[k].bham.a; BHAMB(nbfp,atnr,i,j) = idef->iparams[k].bham.b; BHAMC(nbfp,atnr,i,j) = idef->iparams[k].bham.c; } } } else { snew(nbfp,2*atnr*atnr); for(i=k=0; (i<atnr); i++) { for(j=0; (j<atnr); j++,k++) { C6(nbfp,atnr,i,j) = idef->iparams[k].lj.c6; C12(nbfp,atnr,i,j) = idef->iparams[k].lj.c12; } } } return nbfp; }
static void upd_nbfpbu(FILE *log,real *nbfp,int atnr, real fa[],real fb[],real fc[]) { int n,m,k; /* Update the nonbonded force parameters */ for(k=n=0; (n<atnr); n++) { for(m=0; (m<atnr); m++,k++) { BHAMA(nbfp,atnr,n,m) *= fa[k]; BHAMB(nbfp,atnr,n,m) *= fb[k]; BHAMC(nbfp,atnr,n,m) *= fc[k]; } } }
static void pr_nbfp(FILE *fp,real *nbfp,bool bBHAM,int atnr) { int i,j; if(fp) { for(i=0; (i<atnr); i++) { for(j=0; (j<atnr); j++) { fprintf(fp,"%2d - %2d",i,j); if (bBHAM) fprintf(fp," a=%10g, b=%10g, c=%10g\n",BHAMA(nbfp,atnr,i,j), BHAMB(nbfp,atnr,i,j),BHAMC(nbfp,atnr,i,j)); else fprintf(fp," c6=%10g, c12=%10g\n",C6(nbfp,atnr,i,j), C12(nbfp,atnr,i,j)); } } } }
static void set_bham_b_max(FILE *log,t_forcerec *fr,t_mdatoms *mdatoms) { int i,j,tpi,tpj,ntypes,natoms,*type; real b,bmin; real *nbfp; bool first; if(log) fprintf(log,"Determining largest Buckingham b parameter for table\n"); nbfp = fr->nbfp; ntypes = fr->ntype; type = mdatoms->typeA; natoms = mdatoms->nr; bmin = -1; first = TRUE; fr->bham_b_max = 0; for(i=0; (i<natoms); i++) { tpi = type[i]; if (tpi >= ntypes) fatal_error(0,"Atomtype[%d] = %d, maximum = %d",i,tpi,ntypes); for(j=0; (j<natoms); j++) { tpj = type[j]; if (tpj >= ntypes) fatal_error(0,"Atomtype[%d] = %d, maximum = %d",j,tpj,ntypes); b = BHAMB(nbfp,ntypes,tpi,tpj); if (b > fr->bham_b_max) fr->bham_b_max = b; if ((b < bmin) || first) { bmin = b; first = FALSE; } } } if(log) fprintf(log,"Buckingham b parameters, min: %g, max: %g\n", bmin,fr->bham_b_max); }
void do_coupling(FILE *log,int nfile,t_filenm fnm[], t_coupl_rec *tcr,real t,int step,real ener[], t_forcerec *fr,t_inputrec *ir,bool bMaster, t_mdatoms *md,t_idef *idef,real mu_aver,int nmols, t_commrec *cr,matrix box,tensor virial, tensor pres,rvec mu_tot, rvec x[],rvec f[],bool bDoIt) { #define enm2Debye 48.0321 #define d2e(x) (x)/enm2Debye #define enm2kjmol(x) (x)*0.0143952 /* = 2.0*4.0*M_PI*EPSILON0 */ static real *f6,*f12,*fa,*fb,*fc,*fq; static bool bFirst = TRUE; int i,j,ati,atj,atnr2,type,ftype; real deviation[eoObsNR],prdev[eoObsNR],epot0,dist,rmsf; real ff6,ff12,ffa,ffb,ffc,ffq,factor,dt,mu_ind; real Epol,Eintern,Virial,muabs,xiH=-1,xiS=-1,xi6,xi12; rvec fmol[2]; bool bTest,bPrint; t_coupl_LJ *tclj; t_coupl_BU *tcbu; t_coupl_Q *tcq; t_coupl_iparams *tip; atnr2 = idef->atnr * idef->atnr; if (bFirst) { if (PAR(cr)) fprintf(log,"GCT: this is parallel\n"); else fprintf(log,"GCT: this is not parallel\n"); fflush(log); snew(f6, atnr2); snew(f12,atnr2); snew(fa, atnr2); snew(fb, atnr2); snew(fc, atnr2); snew(fq, idef->atnr); if (tcr->bVirial) { int nrdf = 0; real TTT = 0; real Vol = det(box); for(i=0; (i<ir->opts.ngtc); i++) { nrdf += ir->opts.nrdf[i]; TTT += ir->opts.nrdf[i]*ir->opts.ref_t[i]; } TTT /= nrdf; /* Calculate reference virial from reference temperature and pressure */ tcr->ref_value[eoVir] = 0.5*BOLTZ*nrdf*TTT - (3.0/2.0)* Vol*tcr->ref_value[eoPres]; fprintf(log,"GCT: TTT = %g, nrdf = %d, vir0 = %g, Vol = %g\n", TTT,nrdf,tcr->ref_value[eoVir],Vol); fflush(log); } bFirst = FALSE; } bPrint = MASTER(cr) && do_per_step(step,ir->nstlog); dt = ir->delta_t; /* Initiate coupling to the reference pressure and temperature to start * coupling slowly. */ if (step == 0) { for(i=0; (i<eoObsNR); i++) tcr->av_value[i] = tcr->ref_value[i]; if ((tcr->ref_value[eoDipole]) != 0.0) { mu_ind = mu_aver - d2e(tcr->ref_value[eoDipole]); /* in e nm */ Epol = mu_ind*mu_ind/(enm2kjmol(tcr->ref_value[eoPolarizability])); tcr->av_value[eoEpot] -= Epol; fprintf(log,"GCT: mu_aver = %g(D), mu_ind = %g(D), Epol = %g (kJ/mol)\n", mu_aver*enm2Debye,mu_ind*enm2Debye,Epol); } } /* We want to optimize the LJ params, usually to the Vaporization energy * therefore we only count intermolecular degrees of freedom. * Note that this is now optional. switch UseEinter to yes in your gct file * if you want this. */ dist = calc_dist(log,x); muabs = norm(mu_tot); Eintern = Ecouple(tcr,ener)/nmols; Virial = virial[XX][XX]+virial[YY][YY]+virial[ZZ][ZZ]; /*calc_force(md->nr,f,fmol);*/ clear_rvec(fmol[0]); /* Use a memory of tcr->nmemory steps, so we actually couple to the * average observable over the last tcr->nmemory steps. This may help * in avoiding local minima in parameter space. */ set_act_value(tcr,eoPres, ener[F_PRES],step); set_act_value(tcr,eoEpot, Eintern, step); set_act_value(tcr,eoVir, Virial, step); set_act_value(tcr,eoDist, dist, step); set_act_value(tcr,eoMu, muabs, step); set_act_value(tcr,eoFx, fmol[0][XX], step); set_act_value(tcr,eoFy, fmol[0][YY], step); set_act_value(tcr,eoFz, fmol[0][ZZ], step); set_act_value(tcr,eoPx, pres[XX][XX],step); set_act_value(tcr,eoPy, pres[YY][YY],step); set_act_value(tcr,eoPz, pres[ZZ][ZZ],step); epot0 = tcr->ref_value[eoEpot]; /* If dipole != 0.0 assume we want to use polarization corrected coupling */ if ((tcr->ref_value[eoDipole]) != 0.0) { mu_ind = mu_aver - d2e(tcr->ref_value[eoDipole]); /* in e nm */ Epol = mu_ind*mu_ind/(enm2kjmol(tcr->ref_value[eoPolarizability])); epot0 -= Epol; if (debug) { fprintf(debug,"mu_ind: %g (%g D) mu_aver: %g (%g D)\n", mu_ind,mu_ind*enm2Debye,mu_aver,mu_aver*enm2Debye); fprintf(debug,"Eref %g Epol %g Erunav %g Eact %g\n", tcr->ref_value[eoEpot],Epol,tcr->av_value[eoEpot], tcr->act_value[eoEpot]); } } if (bPrint) { pr_ff(tcr,t,idef,cr,nfile,fnm); } /* Calculate the deviation of average value from the target value */ for(i=0; (i<eoObsNR); i++) { deviation[i] = calc_deviation(tcr->av_value[i],tcr->act_value[i], tcr->ref_value[i]); prdev[i] = tcr->ref_value[i] - tcr->act_value[i]; } deviation[eoEpot] = calc_deviation(tcr->av_value[eoEpot],tcr->act_value[eoEpot], epot0); prdev[eoEpot] = epot0 - tcr->act_value[eoEpot]; if (bPrint) pr_dev(tcr,t,prdev,cr,nfile,fnm); /* First set all factors to 1 */ for(i=0; (i<atnr2); i++) { f6[i] = f12[i] = fa[i] = fb[i] = fc[i] = 1.0; } for(i=0; (i<idef->atnr); i++) fq[i] = 1.0; /* Now compute the actual coupling compononents */ if (!fr->bBHAM) { if (bDoIt) { for(i=0; (i<tcr->nLJ); i++) { tclj=&(tcr->tcLJ[i]); ati=tclj->at_i; atj=tclj->at_j; ff6 = ff12 = 1.0; if (tclj->eObs == eoForce) { gmx_fatal(FARGS,"Hack code for this to work again "); if (debug) fprintf(debug,"Have computed derivatives: xiH = %g, xiS = %g\n",xiH,xiS); if (ati == 1) { /* Hydrogen */ ff12 += xiH; } else if (ati == 2) { /* Shell */ ff12 += xiS; } else gmx_fatal(FARGS,"No H, no Shell, edit code at %s, line %d\n", __FILE__,__LINE__); if (ff6 > 0) set_factor_matrix(idef->atnr,f6, sqrt(ff6), ati,atj); if (ff12 > 0) set_factor_matrix(idef->atnr,f12,sqrt(ff12),ati,atj); } else { if (debug) fprintf(debug,"tcr->tcLJ[%d].xi_6 = %g, xi_12 = %g deviation = %g\n",i, tclj->xi_6,tclj->xi_12,deviation[tclj->eObs]); factor=deviation[tclj->eObs]; upd_f_value(log,idef->atnr,tclj->xi_6, dt,factor,f6, ati,atj); upd_f_value(log,idef->atnr,tclj->xi_12,dt,factor,f12,ati,atj); } } } if (PAR(cr)) { gprod(cr,atnr2,f6); gprod(cr,atnr2,f12); #ifdef DEBUGGCT dump_fm(log,idef->atnr,f6,"f6"); dump_fm(log,idef->atnr,f12,"f12"); #endif } upd_nbfplj(log,fr->nbfp,idef->atnr,f6,f12,tcr->combrule); /* Copy for printing */ for(i=0; (i<tcr->nLJ); i++) { tclj=&(tcr->tcLJ[i]); ati = tclj->at_i; atj = tclj->at_j; if (atj == -1) atj = ati; tclj->c6 = C6(fr->nbfp,fr->ntype,ati,atj); tclj->c12 = C12(fr->nbfp,fr->ntype,ati,atj); } } else { if (bDoIt) { for(i=0; (i<tcr->nBU); i++) { tcbu = &(tcr->tcBU[i]); factor = deviation[tcbu->eObs]; ati = tcbu->at_i; atj = tcbu->at_j; upd_f_value(log,idef->atnr,tcbu->xi_a,dt,factor,fa,ati,atj); upd_f_value(log,idef->atnr,tcbu->xi_b,dt,factor,fb,ati,atj); upd_f_value(log,idef->atnr,tcbu->xi_c,dt,factor,fc,ati,atj); } } if (PAR(cr)) { gprod(cr,atnr2,fa); gprod(cr,atnr2,fb); gprod(cr,atnr2,fc); } upd_nbfpbu(log,fr->nbfp,idef->atnr,fa,fb,fc); /* Copy for printing */ for(i=0; (i<tcr->nBU); i++) { tcbu=&(tcr->tcBU[i]); ati = tcbu->at_i; atj = tcbu->at_j; if (atj == -1) atj = ati; tcbu->a = BHAMA(fr->nbfp,fr->ntype,ati,atj); tcbu->b = BHAMB(fr->nbfp,fr->ntype,ati,atj); tcbu->c = BHAMC(fr->nbfp,fr->ntype,ati,atj); if (debug) fprintf(debug,"buck (type=%d) = %e, %e, %e\n", tcbu->at_i,tcbu->a,tcbu->b,tcbu->c); } } if (bDoIt) { for(i=0; (i<tcr->nQ); i++) { tcq=&(tcr->tcQ[i]); if (tcq->xi_Q) ffq = 1.0 + (dt/tcq->xi_Q) * deviation[tcq->eObs]; else ffq = 1.0; fq[tcq->at_i] *= ffq; } } if (PAR(cr)) gprod(cr,idef->atnr,fq); for(j=0; (j<md->nr); j++) { md->chargeA[j] *= fq[md->typeA[j]]; } for(i=0; (i<tcr->nQ); i++) { tcq=&(tcr->tcQ[i]); for(j=0; (j<md->nr); j++) { if (md->typeA[j] == tcq->at_i) { tcq->Q = md->chargeA[j]; break; } } if (j == md->nr) gmx_fatal(FARGS,"Coupling type %d not found",tcq->at_i); } for(i=0; (i<tcr->nIP); i++) { tip = &(tcr->tIP[i]); type = tip->type; ftype = idef->functype[type]; factor = dt*deviation[tip->eObs]; switch(ftype) { case F_BONDS: if (tip->xi.harmonic.krA) idef->iparams[type].harmonic.krA *= (1+factor/tip->xi.harmonic.krA); if (tip->xi.harmonic.rA) idef->iparams[type].harmonic.rA *= (1+factor/tip->xi.harmonic.rA); break; default: break; } tip->iprint=idef->iparams[type]; } }
static void update_ff(t_forcerec *fr,int nparm,t_range range[],int param_val[]) { static double *sigma=NULL,*eps=NULL,*c6=NULL,*cn=NULL,*bhama=NULL,*bhamb=NULL,*bhamc=NULL; real val,*nbfp; int i,j,atnr; atnr = fr->ntype; nbfp = fr->nbfp; if (fr->bBHAM) { if (bhama == NULL) { snew(bhama,atnr); snew(bhamb,atnr); snew(bhamc,atnr); } } else { if (sigma == NULL) { snew(sigma,atnr); snew(eps,atnr); snew(c6,atnr); snew(cn,atnr); } } /* Get current values for everything */ for(i=0; (i<nparm); i++) { if (ga) val = range[i].rval; else val = value_range(&range[i],param_val[i]); if(debug) fprintf(debug,"val = %g\n",val); switch (range[i].ptype) { case eseSIGMA: sigma[range[i].atype] = val; break; case eseEPSILON: eps[range[i].atype] = val; break; case eseBHAMA: bhama[range[i].atype] = val; break; case eseBHAMB: bhamb[range[i].atype] = val; break; case eseBHAMC: bhamc[range[i].atype] = val; break; case eseCELLX: scale[XX] = val; break; case eseCELLY: scale[YY] = val; break; case eseCELLZ: scale[ZZ] = val; break; default: gmx_fatal(FARGS,"Unknown ptype"); } } if (fr->bBHAM) { for(i=0; (i<atnr); i++) { for(j=0; (j<=i); j++) { BHAMA(nbfp,atnr,i,j) = BHAMA(nbfp,atnr,j,i) = sqrt(bhama[i]*bhama[j]); BHAMB(nbfp,atnr,i,j) = BHAMB(nbfp,atnr,j,i) = sqrt(bhamb[i]*bhamb[j]); BHAMC(nbfp,atnr,i,j) = BHAMC(nbfp,atnr,j,i) = sqrt(bhamc[i]*bhamc[j]); } } } else { /* Now build a new matrix */ for(i=0; (i<atnr); i++) { c6[i] = 4*eps[i]*pow(sigma[i],6.0); cn[i] = 4*eps[i]*pow(sigma[i],ff.npow); } for(i=0; (i<atnr); i++) { for(j=0; (j<=i); j++) { C6(nbfp,atnr,i,j) = C6(nbfp,atnr,j,i) = sqrt(c6[i]*c6[j]); C12(nbfp,atnr,i,j) = C12(nbfp,atnr,j,i) = sqrt(cn[i]*cn[j]); } } } if (debug) { if (!fr->bBHAM) for(i=0; (i<atnr); i++) fprintf(debug,"atnr = %2d sigma = %8.4f eps = %8.4f\n",i,sigma[i],eps[i]); for(i=0; (i<atnr); i++) { for(j=0; (j<atnr); j++) { if (fr->bBHAM) fprintf(debug,"i: %2d j: %2d A: %10.5e B: %10.5e C: %10.5e\n",i,j, BHAMA(nbfp,atnr,i,j),BHAMB(nbfp,atnr,i,j),BHAMC(nbfp,atnr,i,j)); else fprintf(debug,"i: %2d j: %2d c6: %10.5e cn: %10.5e\n",i,j, C6(nbfp,atnr,i,j),C12(nbfp,atnr,i,j)); } } } }
static void check_solvent(FILE *fp,t_topology *top,t_forcerec *fr, t_mdatoms *md,t_nsborder *nsb) { /* This routine finds out whether a charge group can be used as * solvent in the innerloops. The routine should be called once * at the beginning of the MD program. */ t_block *cgs,*excl,*mols; atom_id *cgid; int i,j,m,j0,j1,nj,k,aj,ak,tjA,tjB,nl_m,nl_n,nl_o; int warncount; bool bOneCG; bool *bAllExcl,bAE,bOrder; bool *bHaveLJ,*bHaveCoul; cgs = &(top->blocks[ebCGS]); excl = &(top->atoms.excl); mols = &(top->blocks[ebMOLS]); if (fp) fprintf(fp,"Going to determine what solvent types we have.\n"); snew(fr->solvent_type,cgs->nr+1); snew(fr->mno_index,(cgs->nr+1)*3); /* Generate charge group number for all atoms */ cgid = make_invblock(cgs,cgs->nra); warncount=0; /* Loop over molecules */ if (fp) fprintf(fp,"There are %d molecules, %d charge groups and %d atoms\n", mols->nr,cgs->nr,cgs->nra); for(i=0; (i<mols->nr); i++) { /* Set boolean that determines whether the molecules consists of one CG */ bOneCG = TRUE; /* Set some counters */ j0 = mols->index[i]; j1 = mols->index[i+1]; nj = j1-j0; for(j=j0+1; (j<j1); j++) { bOneCG = bOneCG && (cgid[mols->a[j]] == cgid[mols->a[j-1]]); } if (fr->bSolvOpt && bOneCG && nj>1) { /* Check whether everything is excluded */ snew(bAllExcl,nj); bAE = TRUE; /* Loop over all atoms in molecule */ for(j=j0; (j<j1) && bAE; j++) { /* Set a flag for each atom in the molecule that determines whether * it is excluded or not */ for(k=0; (k<nj); k++) bAllExcl[k] = FALSE; /* Now check all the exclusions of this atom */ for(k=excl->index[j]; (k<excl->index[j+1]); k++) { ak = excl->a[k]; /* Consistency and range check */ if ((ak < j0) || (ak >= j1)) fatal_error(0,"Exclusion outside molecule? ak = %d, j0 = %d, j1 = 5d, mol is %d",ak,j0,j1,i); bAllExcl[ak-j0] = TRUE; } /* Now sum up the booleans */ for(k=0; (k<nj); k++) bAE = bAE && bAllExcl[k]; } if (bAE) { snew(bHaveCoul,nj); snew(bHaveLJ,nj); for(j=j0; (j<j1); j++) { /* Check for coulomb */ aj = mols->a[j]; bHaveCoul[j-j0] = ( (fabs(top->atoms.atom[aj].q ) > GMX_REAL_MIN) || (fabs(top->atoms.atom[aj].qB) > GMX_REAL_MIN)); /* Check for LJ. */ tjA = top->atoms.atom[aj].type; tjB = top->atoms.atom[aj].typeB; bHaveLJ[j-j0] = FALSE; for(k=0; (k<fr->ntype); k++) { if (fr->bBHAM) bHaveLJ[j-j0] = (bHaveLJ[j-j0] || (fabs(BHAMA(fr->nbfp,fr->ntype,tjA,k)) > GMX_REAL_MIN) || (fabs(BHAMB(fr->nbfp,fr->ntype,tjA,k)) > GMX_REAL_MIN) || (fabs(BHAMC(fr->nbfp,fr->ntype,tjA,k)) > GMX_REAL_MIN) || (fabs(BHAMA(fr->nbfp,fr->ntype,tjB,k)) > GMX_REAL_MIN) || (fabs(BHAMB(fr->nbfp,fr->ntype,tjB,k)) > GMX_REAL_MIN) || (fabs(BHAMC(fr->nbfp,fr->ntype,tjB,k)) > GMX_REAL_MIN)); else bHaveLJ[j-j0] = (bHaveLJ[j-j0] || (fabs(C6(fr->nbfp,fr->ntype,tjA,k)) > GMX_REAL_MIN) || (fabs(C12(fr->nbfp,fr->ntype,tjA,k)) > GMX_REAL_MIN) || (fabs(C6(fr->nbfp,fr->ntype,tjB,k)) > GMX_REAL_MIN) || (fabs(C12(fr->nbfp,fr->ntype,tjB,k)) > GMX_REAL_MIN)); } } /* Now we have determined what particles have which interactions * In the case of water-like molecules we only check for the number * of particles and the LJ, not for the Coulomb. Let's just assume * that the water loops are faster than the MNO loops anyway. DvdS */ /* No - there's another problem: To optimize the water * innerloop assumes the charge of the first i atom is constant * qO, and charge on atoms 2/3 is constant qH. /EL */ /* I won't write any altivec versions of the general solvent inner * loops. Thus, when USE_PPC_ALTIVEC is defined it is faster * to use the normal loops instead of the MNO solvent version. /EL */ aj=mols->a[j0]; if((nj==3) && bHaveCoul[0] && bHaveLJ[0] && !bHaveLJ[1] && !bHaveLJ[2] && fabs(top->atoms.atom[aj+1].q - top->atoms.atom[aj+2].q) < GMX_REAL_MIN) fr->solvent_type[cgid[aj]] = esolWATER; else { #ifdef USE_PPC_ALTIVEC fr->solvent_type[cgid[aj]] = esolNO; #else /* Time to compute M & N & O */ for(k=0; (k<nj) && (bHaveLJ[k] && bHaveCoul[k]); k++) ; nl_n = k; for(; (k<nj) && (!bHaveLJ[k] && bHaveCoul[k]); k++) ; nl_o = k; for(; (k<nj) && (bHaveLJ[k] && !bHaveCoul[k]); k++) ; nl_m = k; /* Now check whether we're at the end of the pack */ bOrder = FALSE; for(; (k<nj); k++) bOrder = bOrder || (bHaveLJ[k] || bHaveCoul[k]); if (bOrder) { /* If we have a solvent molecule with LJC everywhere, then * we shouldn't issue a warning. Only if we suspect something * could be better. */ if (nl_n != nj) { warncount++; if(warncount<11 && fp) fprintf(fp,"The order in molecule %d could be optimized" " for better performance\n",i); if(warncount==10 && fp) fprintf(fp,"(More than 10 molecules where the order can be optimized)\n"); } nl_m = nl_n = nl_o = nj; } fr->mno_index[cgid[aj]*3] = nl_m; fr->mno_index[cgid[aj]*3+1] = nl_n; fr->mno_index[cgid[aj]*3+2] = nl_o; fr->solvent_type[cgid[aj]] = esolMNO; #endif /* MNO solvent if not using altivec */ } /* Last check for perturbed atoms */ for(j=j0; (j<j1); j++) if (md->bPerturbed[mols->a[j]]) fr->solvent_type[cgid[mols->a[j0]]] = esolNO; sfree(bHaveLJ); sfree(bHaveCoul); } else { /* Turn off solvent optimization for all cg's in the molecule, * here there is only one. */ fr->solvent_type[cgid[mols->a[j0]]] = esolNO; } sfree(bAllExcl); } else { /* Turn off solvent optimization for all cg's in the molecule */ for(j=mols->index[i]; (j<mols->index[i+1]); j++) { fr->solvent_type[cgid[mols->a[j]]] = esolNO; } } } if (debug) { for(i=0; (i<cgs->nr); i++) fprintf(debug,"MNO: cg = %5d, m = %2d, n = %2d, o = %2d\n", i,fr->mno_index[3*i],fr->mno_index[3*i+1],fr->mno_index[3*i+2]); } /* Now compute the number of solvent molecules, could be merged with code above */ fr->nMNOMol = 0; fr->nWatMol = 0; for(m=0; m<3; m++) fr->nMNOav[m] = 0; for(i=0; i<mols->nr; i++) { j = mols->a[mols->index[i]]; if (j>=START(nsb) && j<START(nsb)+HOMENR(nsb)) { if (fr->solvent_type[cgid[j]] == esolMNO) { fr->nMNOMol++; for(m=0; m<3; m++) fr->nMNOav[m] += fr->mno_index[3*cgid[j]+m]; } else if (fr->solvent_type[cgid[j]] == esolWATER) fr->nWatMol++; } } if (fr->nMNOMol > 0) for(m=0; (m<3); m++) fr->nMNOav[m] /= fr->nMNOMol; sfree(cgid); if (fp) { fprintf(fp,"There are %d optimized solvent molecules on node %d\n", fr->nMNOMol,nsb->nodeid); if (fr->nMNOMol > 0) fprintf(fp," aver. nr. of atoms per molecule: vdwc %.1f coul %.1f vdw %.1f\n", fr->nMNOav[1],fr->nMNOav[2]-fr->nMNOav[1],fr->nMNOav[0]-fr->nMNOav[2]); fprintf(fp,"There are %d optimized water molecules on node %d\n", fr->nWatMol,nsb->nodeid); } }