1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
50 #include "nonbonded.h"
79 gmx_grppairener_t *grppener,
81 gmx_bool bDoLongRange,
89 if (!fr->ns.nblist_initialized)
91 init_neighbor_list(fp, fr, md->homenr);
97 nsearch = search_neighbours(fp,fr,x,box,top,groups,cr,nrnb,md,
98 lambda,dvdlambda,grppener,
99 bFillGrid,bDoLongRange,
102 fprintf(debug,"nsearch = %d\n",nsearch);
104 /* Check whether we have to do dynamic load balancing */
105 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
106 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
107 &(top->idef),opts->ngener);
109 if (fr->ns.dump_nl > 0)
110 dump_nblist(fp,cr,fr,fr->ns.dump_nl);
113 void do_force_lowlevel(FILE *fplog, gmx_large_int_t step,
114 t_forcerec *fr, t_inputrec *ir,
115 t_idef *idef, t_commrec *cr,
116 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
119 rvec x[], history_t *hist,
121 gmx_enerdata_t *enerd,
139 gmx_bool bDoEpot,bSepDVDL,bSB;
143 real Vsr,Vlr,Vcorr=0,vdip,vcharge;
147 gmx_enerdata_t ed_lam;
148 double clam_i,vlam_i;
149 real dvdl_dum[efptNR], dvdlambda[efptNR], lam_i[efptNR];
150 real dvdlsum,dvdl_walls;
153 double t0=0.0,t1,t2,t3; /* time measurement for coarse load balancing */
156 #define PRINT_SEPDVDL(s,v,dvdlambda) if (bSepDVDL) fprintf(fplog,sepdvdlformat,s,v,dvdlambda);
159 set_pbc(&pbc,fr->ePBC,box);
161 /* reset free energy components */
162 for (i=0;i<efptNR;i++)
169 for(i=0; (i<DIM); i++)
171 box_size[i]=box[i][i];
174 bSepDVDL=(fr->bSepDVDL && do_per_step(step,ir->nstlog));
177 /* do QMMM first if requested */
180 enerd->term[F_EQM] = calculate_QMMM(cr,x,f,fr,md);
185 fprintf(fplog,"Step %s: non-bonded V and dVdl for node %d:\n",
186 gmx_step_str(step,buf),cr->nodeid);
189 /* Call the short range functions all in one go. */
192 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
193 #define TAKETIME FALSE
196 MPI_Barrier(cr->mpi_comm_mygroup);
203 /* foreign lambda component for walls */
204 dvdl_walls = do_walls(ir,fr,box,md,x,f,lambda[efptVDW],
205 enerd->grpp.ener[egLJSR],nrnb);
206 PRINT_SEPDVDL("Walls",0.0,dvdl_walls);
207 dvdlambda[efptVDW] += dvdl_walls;
208 enerd->dvdl_lin[efptVDW] += dvdl_walls;
211 /* If doing GB, reset dvda and calculate the Born radii */
212 if (ir->implicit_solvent)
214 /* wallcycle_start(wcycle,ewcGB); */
216 for(i=0;i<born->nr;i++)
223 calc_gb_rad(cr,fr,ir,top,atype,x,&(fr->gblist),born,md,nrnb);
226 /* wallcycle_stop(wcycle, ewcGB); */
231 if (flags & GMX_FORCE_FORCES)
233 donb_flags |= GMX_DONB_FORCES;
236 do_nonbonded(cr,fr,x,f,md,excl,
238 enerd->grpp.ener[egBHAMSR] :
239 enerd->grpp.ener[egLJSR],
240 enerd->grpp.ener[egCOULSR],
241 enerd->grpp.ener[egGB],box_size,nrnb,
242 lambda,dvdlambda,-1,-1,donb_flags);
243 /* If we do foreign lambda and we have soft-core interactions
244 * we have to recalculate the (non-linear) energies contributions.
246 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
248 init_enerdata(mtop->groups.grps[egcENER].nr,fepvals->n_lambda,&ed_lam);
250 for(i=0; i<enerd->n_lambda; i++)
252 for (j=0;j<efptNR;j++)
254 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
256 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
257 do_nonbonded(cr,fr,x,f,md,excl,
259 ed_lam.grpp.ener[egBHAMSR] :
260 ed_lam.grpp.ener[egLJSR],
261 ed_lam.grpp.ener[egCOULSR],
262 enerd->grpp.ener[egGB], box_size,nrnb,
263 lam_i,dvdl_dum,-1,-1,
264 GMX_DONB_FOREIGNLAMBDA);
265 sum_epot(&ir->opts,&ed_lam);
266 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
268 destroy_enerdata(&ed_lam);
272 /* If we are doing GB, calculate bonded forces and apply corrections
273 * to the solvation forces */
274 /* MRS: Eventually, many need to include free energy contribution here! */
275 if (ir->implicit_solvent) {
276 calc_gb_forces(cr,md,born,top,atype,x,f,fr,idef,
277 ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&pbc,graph,enerd);
288 if (fepvals->sc_alpha!=0)
290 enerd->dvdl_nonlin[efptVDW] += dvdlambda[efptVDW];
294 enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
297 if (fepvals->sc_alpha!=0)
299 /* even though coulomb part is linear, we already added it, beacuse we
300 need to go through the vdw calculation anyway */
302 enerd->dvdl_nonlin[efptCOUL] += dvdlambda[efptCOUL];
306 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
312 for(i=0; i<enerd->grpp.nener; i++)
316 enerd->grpp.ener[egBHAMSR][i] :
317 enerd->grpp.ener[egLJSR][i])
318 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
320 dvdlsum = dvdlambda[efptVDW]+dvdlambda[efptCOUL];
321 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",Vsr,dvdlsum);
328 pr_rvecs(debug,0,"fshift after SR",fr->fshift,SHIFTS);
331 /* Shift the coordinates. Must be done before bonded forces and PPPM,
332 * but is also necessary for SHAKE and update, therefore it can NOT
333 * go when no bonded forces have to be evaluated.
336 /* Here sometimes we would not need to shift with NBFonly,
337 * but we do so anyhow for consistency of the returned coordinates.
341 shift_self(graph,box,x);
344 inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
348 inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
351 /* Check whether we need to do bondeds or correct for exclusions */
353 ((flags & GMX_FORCE_BONDED)
354 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
356 /* Since all atoms are in the rectangular or triclinic unit-cell,
357 * only single box vector shifts (2 in x) are required.
359 set_pbc_dd(&pbc,fr->ePBC,cr->dd,TRUE,box);
363 if (flags & GMX_FORCE_BONDED)
365 calc_bonds(fplog,cr->ms,
366 idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
367 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
368 fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
370 /* Check if we have to determine energy differences
371 * at foreign lambda's.
373 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
374 idef->ilsort != ilsortNO_FE)
376 if (idef->ilsort != ilsortFE_SORTED)
378 gmx_incons("The bonded interactions are not sorted for free energy");
380 init_enerdata(mtop->groups.grps[egcENER].nr,fepvals->n_lambda,&ed_lam);
382 for(i=0; i<enerd->n_lambda; i++)
384 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
385 for (j=0;j<efptNR;j++)
387 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
389 calc_bonds_lambda(fplog,idef,x,fr,&pbc,graph,&ed_lam,nrnb,lam_i,md,
390 fcd,DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
391 sum_epot(&ir->opts,&ed_lam);
392 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
394 destroy_enerdata(&ed_lam);
402 if (EEL_FULL(fr->eeltype))
404 bSB = (ir->nwall == 2);
408 svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
409 box_size[ZZ] *= ir->wall_ewald_zfac;
412 clear_mat(fr->vir_el_recip);
418 dvdlambda[efptCOUL] = 0;
419 Vcorr = ewald_LRcorrection(fplog,md->start,md->start+md->homenr,
422 md->nChargePerturbed ? md->chargeB : NULL,
423 excl,x,bSB ? boxs : box,mu_tot,
426 lambda[efptCOUL],&dvdlambda[efptCOUL],&vdip,&vcharge);
427 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdlambda);
428 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
432 if (ir->ewald_geometry != eewg3D || ir->epsilon_surface != 0)
434 gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
436 /* The TPI molecule does not have exclusions with the rest
437 * of the system and no intra-molecular PME grid contributions
438 * will be calculated in gmx_pme_calc_energy.
444 dvdlambda[efptCOUL] = 0;
451 case eelPMEUSERSWITCH:
453 if (cr->duty & DUTY_PME)
455 assert(fr->n_tpi >= 0);
456 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
458 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
459 if (flags & GMX_FORCE_FORCES)
461 pme_flags |= GMX_PME_CALC_F;
463 if (flags & GMX_FORCE_VIRIAL)
465 pme_flags |= GMX_PME_CALC_ENER_VIR;
469 /* We don't calculate f, but we do want the potential */
470 pme_flags |= GMX_PME_CALC_POT;
472 wallcycle_start(wcycle,ewcPMEMESH);
473 status = gmx_pme_do(fr->pmedata,
474 md->start,md->homenr - fr->n_tpi,
476 md->chargeA,md->chargeB,
478 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
479 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
481 fr->vir_el_recip,fr->ewaldcoeff,
482 &Vlr,lambda[efptCOUL],&dvdlambda[efptCOUL],
484 *cycles_pme = wallcycle_stop(wcycle,ewcPMEMESH);
486 /* We should try to do as little computation after
487 * this as possible, because parallel PME synchronizes
488 * the nodes, so we want all load imbalance of the rest
489 * of the force calculation to be before the PME call.
490 * DD load balancing is done on the whole time of
491 * the force call (without PME).
496 /* Determine the PME grid energy of the test molecule
497 * with the PME grid potential of the other charges.
499 gmx_pme_calc_energy(fr->pmedata,fr->n_tpi,
500 x + md->homenr - fr->n_tpi,
501 md->chargeA + md->homenr - fr->n_tpi,
504 PRINT_SEPDVDL("PME mesh",Vlr,dvdlambda[efptCOUL]);
508 /* Energies and virial are obtained later from the PME nodes */
509 /* but values have to be zeroed out here */
514 Vlr = do_ewald(fplog,FALSE,ir,x,fr->f_novirsum,
515 md->chargeA,md->chargeB,
516 box_size,cr,md->homenr,
517 fr->vir_el_recip,fr->ewaldcoeff,
518 lambda[efptCOUL],&dvdlambda[efptCOUL],fr->ewald_table);
519 PRINT_SEPDVDL("Ewald long-range",Vlr,dvdlambda[efptCOUL]);
523 gmx_fatal(FARGS,"No such electrostatics method implemented %s",
524 eel_names[fr->eeltype]);
528 gmx_fatal(FARGS,"Error %d in long range electrostatics routine %s",
529 status,EELTYPE(fr->eeltype));
531 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
532 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
535 fprintf(debug,"Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
536 Vlr,Vcorr,enerd->term[F_COUL_RECIP]);
537 pr_rvecs(debug,0,"vir_el_recip after corr",fr->vir_el_recip,DIM);
538 pr_rvecs(debug,0,"fshift after LR Corrections",fr->fshift,SHIFTS);
543 if (EEL_RF(fr->eeltype))
545 dvdlambda[efptCOUL] = 0;
547 if (fr->eeltype != eelRF_NEC)
549 enerd->term[F_RF_EXCL] =
550 RF_excl_correction(fplog,fr,graph,md,excl,x,f,
551 fr->fshift,&pbc,lambda[efptCOUL],&dvdlambda[efptCOUL]);
554 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
555 PRINT_SEPDVDL("RF exclusion correction",
556 enerd->term[F_RF_EXCL],dvdlambda[efptCOUL]);
564 print_nrnb(debug,nrnb);
572 MPI_Barrier(cr->mpi_comm_mygroup);
575 if (fr->timesteps == 11)
577 fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
578 cr->nodeid, gmx_step_str(fr->timesteps,buf),
579 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
580 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
588 pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
593 void init_enerdata(int ngener,int n_lambda,gmx_enerdata_t *enerd)
597 for(i=0; i<F_NRE; i++)
603 for(i=0; i<efptNR; i++) {
604 enerd->dvdl_lin[i] = 0;
605 enerd->dvdl_nonlin[i] = 0;
611 fprintf(debug,"Creating %d sized group matrix for energies\n",n2);
613 enerd->grpp.nener = n2;
614 for(i=0; (i<egNR); i++)
616 snew(enerd->grpp.ener[i],n2);
621 enerd->n_lambda = 1 + n_lambda;
622 snew(enerd->enerpart_lambda,enerd->n_lambda);
630 void destroy_enerdata(gmx_enerdata_t *enerd)
634 for(i=0; (i<egNR); i++)
636 sfree(enerd->grpp.ener[i]);
641 sfree(enerd->enerpart_lambda);
645 static real sum_v(int n,real v[])
657 void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd)
659 gmx_grppairener_t *grpp;
666 /* Accumulate energies */
667 epot[F_COUL_SR] = sum_v(grpp->nener,grpp->ener[egCOULSR]);
668 epot[F_LJ] = sum_v(grpp->nener,grpp->ener[egLJSR]);
669 epot[F_LJ14] = sum_v(grpp->nener,grpp->ener[egLJ14]);
670 epot[F_COUL14] = sum_v(grpp->nener,grpp->ener[egCOUL14]);
671 epot[F_COUL_LR] = sum_v(grpp->nener,grpp->ener[egCOULLR]);
672 epot[F_LJ_LR] = sum_v(grpp->nener,grpp->ener[egLJLR]);
673 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
674 epot[F_GBPOL] += sum_v(grpp->nener,grpp->ener[egGB]);
676 /* lattice part of LR doesnt belong to any group
677 * and has been added earlier
679 epot[F_BHAM] = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
680 epot[F_BHAM_LR] = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
683 for(i=0; (i<F_EPOT); i++)
685 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
687 epot[F_EPOT] += epot[i];
692 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
697 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
698 enerd->term[F_DVDL] = 0.0;
699 for (i=0;i<efptNR;i++)
701 if (fepvals->separate_dvdl[i])
703 /* could this be done more readably/compactly? */
712 index = F_DVDL_BONDED;
714 case (efptRESTRAINT):
715 index = F_DVDL_RESTRAINT;
724 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
727 fprintf(debug,"dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
728 efpt_names[i],i,enerd->term[index],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
733 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
736 fprintf(debug,"dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
737 efpt_names[0],i,enerd->term[F_DVDL],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
742 /* Notes on the foreign lambda free energy difference evaluation:
743 * Adding the potential and ekin terms that depend linearly on lambda
744 * as delta lam * dvdl to the energy differences is exact.
745 * For the constraints this is not exact, but we have no other option
746 * without literally changing the lengths and reevaluating the energies at each step.
747 * (try to remedy this post 4.6 - MRS)
748 * For the non-bonded LR term we assume that the soft-core (if present)
749 * no longer affects the energy beyond the short-range cut-off,
750 * which is a very good approximation (except for exotic settings).
751 * (investigate how to overcome this post 4.6 - MRS)
754 for(i=0; i<fepvals->n_lambda; i++)
755 { /* note we are iterating over fepvals here!
756 For the current lam, dlam = 0 automatically,
757 so we don't need to add anything to the
758 enerd->enerpart_lambda[0] */
760 /* we don't need to worry about dvdl contributions to the current lambda, because
761 it's automatically zero */
763 /* first kinetic energy term */
764 dlam = (fepvals->all_lambda[efptMASS][i] - lambda[efptMASS]);
766 enerd->enerpart_lambda[i+1] += enerd->term[F_DKDL]*dlam;
768 for (j=0;j<efptNR;j++)
770 if (j==efptMASS) {continue;} /* no other mass term to worry about */
772 dlam = (fepvals->all_lambda[j][i]-lambda[j]);
773 enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
776 fprintf(debug,"enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
777 fepvals->all_lambda[j][i],efpt_names[j],
778 (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
779 dlam,enerd->dvdl_lin[j]);
785 void reset_enerdata(t_grpopts *opts,
786 t_forcerec *fr,gmx_bool bNS,
787 gmx_enerdata_t *enerd,
793 /* First reset all energy components, except for the long range terms
794 * on the master at non neighbor search steps, since the long range
795 * terms have already been summed at the last neighbor search step.
797 bKeepLR = (fr->bTwinRange && !bNS);
798 for(i=0; (i<egNR); i++) {
799 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
800 for(j=0; (j<enerd->grpp.nener); j++)
801 enerd->grpp.ener[i][j] = 0.0;
804 for (i=0;i<efptNR;i++)
806 enerd->dvdl_lin[i] = 0.0;
807 enerd->dvdl_nonlin[i] = 0.0;
810 /* Normal potential energy components */
811 for(i=0; (i<=F_EPOT); i++) {
812 enerd->term[i] = 0.0;
814 /* Initialize the dVdlambda term with the long range contribution */
815 /* Initialize the dvdl term with the long range contribution */
816 enerd->term[F_DVDL] = 0.0;
817 enerd->term[F_DVDL_COUL] = 0.0;
818 enerd->term[F_DVDL_VDW] = 0.0;
819 enerd->term[F_DVDL_BONDED] = 0.0;
820 enerd->term[F_DVDL_RESTRAINT] = 0.0;
821 enerd->term[F_DKDL] = 0.0;
822 if (enerd->n_lambda > 0)
824 for(i=0; i<enerd->n_lambda; i++)
826 enerd->enerpart_lambda[i] = 0.0;