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
49 #include "nonbonded.h"
78 gmx_grppairener_t *grppener,
80 gmx_bool bDoLongRange,
88 if (!fr->ns.nblist_initialized)
90 init_neighbor_list(fp, fr, md->homenr);
96 nsearch = search_neighbours(fp,fr,x,box,top,groups,cr,nrnb,md,
97 lambda,dvdlambda,grppener,
98 bFillGrid,bDoLongRange,
101 fprintf(debug,"nsearch = %d\n",nsearch);
103 /* Check whether we have to do dynamic load balancing */
104 /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
105 count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
106 &(top->idef),opts->ngener);
108 if (fr->ns.dump_nl > 0)
109 dump_nblist(fp,cr,fr,fr->ns.dump_nl);
112 void do_force_lowlevel(FILE *fplog, gmx_large_int_t step,
113 t_forcerec *fr, t_inputrec *ir,
114 t_idef *idef, t_commrec *cr,
115 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
118 rvec x[], history_t *hist,
120 gmx_enerdata_t *enerd,
137 gmx_bool bDoEpot,bSepDVDL,bSB;
141 real dvdlambda,Vsr,Vlr,Vcorr=0,vdip,vcharge;
145 gmx_enerdata_t ed_lam;
150 double t0=0.0,t1,t2,t3; /* time measurement for coarse load balancing */
153 #define PRINT_SEPDVDL(s,v,dvdl) if (bSepDVDL) fprintf(fplog,sepdvdlformat,s,v,dvdl);
156 set_pbc(&pbc,fr->ePBC,box);
159 for(i=0; (i<DIM); i++)
161 box_size[i]=box[i][i];
164 bSepDVDL=(fr->bSepDVDL && do_per_step(step,ir->nstlog));
167 /* do QMMM first if requested */
170 enerd->term[F_EQM] = calculate_QMMM(cr,x,f,fr,md);
175 fprintf(fplog,"Step %s: non-bonded V and dVdl for node %d:\n",
176 gmx_step_str(step,buf),cr->nodeid);
179 /* Call the short range functions all in one go. */
184 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
185 #define TAKETIME FALSE
188 MPI_Barrier(cr->mpi_comm_mygroup);
195 dvdlambda = do_walls(ir,fr,box,md,x,f,lambda,
196 enerd->grpp.ener[egLJSR],nrnb);
197 PRINT_SEPDVDL("Walls",0.0,dvdlambda);
198 enerd->dvdl_lin += dvdlambda;
201 /* If doing GB, reset dvda and calculate the Born radii */
202 if (ir->implicit_solvent)
204 /* wallcycle_start(wcycle,ewcGB); */
206 for(i=0;i<born->nr;i++)
213 calc_gb_rad(cr,fr,ir,top,atype,x,&(fr->gblist),born,md,nrnb);
216 /* wallcycle_stop(wcycle, ewcGB); */
221 if (flags & GMX_FORCE_FORCES)
223 donb_flags |= GMX_DONB_FORCES;
225 do_nonbonded(cr,fr,x,f,md,excl,
227 enerd->grpp.ener[egBHAMSR] :
228 enerd->grpp.ener[egLJSR],
229 enerd->grpp.ener[egCOULSR],
230 enerd->grpp.ener[egGB],box_size,nrnb,
231 lambda,&dvdlambda,-1,-1,donb_flags);
232 /* If we do foreign lambda and we have soft-core interactions
233 * we have to recalculate the (non-linear) energies contributions.
235 if (ir->n_flambda > 0 && (flags & GMX_FORCE_DHDL) && ir->sc_alpha != 0)
237 init_enerdata(mtop->groups.grps[egcENER].nr,ir->n_flambda,&ed_lam);
239 for(i=0; i<enerd->n_lambda; i++)
241 lam_i = (i==0 ? lambda : ir->flambda[i-1]);
243 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
244 do_nonbonded(cr,fr,x,f,md,excl,
246 ed_lam.grpp.ener[egBHAMSR] :
247 ed_lam.grpp.ener[egLJSR],
248 ed_lam.grpp.ener[egCOULSR],
249 enerd->grpp.ener[egGB], box_size,nrnb,
250 lam_i,&dvdl_dum,-1,-1,
251 GMX_DONB_FOREIGNLAMBDA);
252 sum_epot(&ir->opts,&ed_lam);
253 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
255 destroy_enerdata(&ed_lam);
259 /* If we are doing GB, calculate bonded forces and apply corrections
260 * to the solvation forces */
261 if (ir->implicit_solvent) {
262 calc_gb_forces(cr,md,born,top,atype,x,f,fr,idef,
263 ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&pbc,graph,enerd);
274 if (ir->sc_alpha != 0)
276 enerd->dvdl_nonlin += dvdlambda;
280 enerd->dvdl_lin += dvdlambda;
285 for(i=0; i<enerd->grpp.nener; i++)
289 enerd->grpp.ener[egBHAMSR][i] :
290 enerd->grpp.ener[egLJSR][i])
291 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
294 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",Vsr,dvdlambda);
299 pr_rvecs(debug,0,"fshift after SR",fr->fshift,SHIFTS);
302 /* Shift the coordinates. Must be done before bonded forces and PPPM,
303 * but is also necessary for SHAKE and update, therefore it can NOT
304 * go when no bonded forces have to be evaluated.
307 /* Here sometimes we would not need to shift with NBFonly,
308 * but we do so anyhow for consistency of the returned coordinates.
312 shift_self(graph,box,x);
315 inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
319 inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
322 /* Check whether we need to do bondeds or correct for exclusions */
324 ((flags & GMX_FORCE_BONDED)
325 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
327 /* Since all atoms are in the rectangular or triclinic unit-cell,
328 * only single box vector shifts (2 in x) are required.
330 set_pbc_dd(&pbc,fr->ePBC,cr->dd,TRUE,box);
334 if (flags & GMX_FORCE_BONDED)
336 calc_bonds(fplog,cr->ms,
337 idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
338 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
339 fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
341 /* Check if we have to determine energy differences
342 * at foreign lambda's.
344 if (ir->n_flambda > 0 && (flags & GMX_FORCE_DHDL) &&
345 idef->ilsort != ilsortNO_FE)
347 if (idef->ilsort != ilsortFE_SORTED)
349 gmx_incons("The bonded interactions are not sorted for free energy");
351 init_enerdata(mtop->groups.grps[egcENER].nr,ir->n_flambda,&ed_lam);
353 for(i=0; i<enerd->n_lambda; i++)
355 lam_i = (i==0 ? lambda : ir->flambda[i-1]);
357 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
358 calc_bonds_lambda(fplog,
359 idef,x,fr,&pbc,graph,&ed_lam,nrnb,lam_i,md,
361 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
362 sum_epot(&ir->opts,&ed_lam);
363 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
365 destroy_enerdata(&ed_lam);
373 if (EEL_FULL(fr->eeltype))
375 bSB = (ir->nwall == 2);
379 svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
380 box_size[ZZ] *= ir->wall_ewald_zfac;
383 clear_mat(fr->vir_el_recip);
390 Vcorr = ewald_LRcorrection(fplog,md->start,md->start+md->homenr,
393 md->nChargePerturbed ? md->chargeB : NULL,
394 excl,x,bSB ? boxs : box,mu_tot,
397 lambda,&dvdlambda,&vdip,&vcharge);
398 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdlambda);
399 enerd->dvdl_lin += dvdlambda;
403 if (ir->ewald_geometry != eewg3D || ir->epsilon_surface != 0)
405 gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
407 /* The TPI molecule does not have exclusions with the rest
408 * of the system and no intra-molecular PME grid contributions
409 * will be calculated in gmx_pme_calc_energy.
422 case eelPMEUSERSWITCH:
424 if (cr->duty & DUTY_PME)
426 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
428 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
429 if (flags & GMX_FORCE_FORCES)
431 pme_flags |= GMX_PME_CALC_F;
433 if (flags & GMX_FORCE_VIRIAL)
435 pme_flags |= GMX_PME_CALC_ENER_VIR;
439 /* We don't calculate f, but we do want the potential */
440 pme_flags |= GMX_PME_CALC_POT;
442 wallcycle_start(wcycle,ewcPMEMESH);
443 status = gmx_pme_do(fr->pmedata,
444 md->start,md->homenr - fr->n_tpi,
446 md->chargeA,md->chargeB,
448 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
449 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
451 fr->vir_el_recip,fr->ewaldcoeff,
452 &Vlr,lambda,&dvdlambda,
454 *cycles_pme = wallcycle_stop(wcycle,ewcPMEMESH);
456 /* We should try to do as little computation after
457 * this as possible, because parallel PME synchronizes
458 * the nodes, so we want all load imbalance of the rest
459 * of the force calculation to be before the PME call.
460 * DD load balancing is done on the whole time of
461 * the force call (without PME).
466 /* Determine the PME grid energy of the test molecule
467 * with the PME grid potential of the other charges.
469 gmx_pme_calc_energy(fr->pmedata,fr->n_tpi,
470 x + md->homenr - fr->n_tpi,
471 md->chargeA + md->homenr - fr->n_tpi,
474 PRINT_SEPDVDL("PME mesh",Vlr,dvdlambda);
478 /* Energies and virial are obtained later from the PME nodes */
479 /* but values have to be zeroed out here */
484 Vlr = do_ewald(fplog,FALSE,ir,x,fr->f_novirsum,
485 md->chargeA,md->chargeB,
486 box_size,cr,md->homenr,
487 fr->vir_el_recip,fr->ewaldcoeff,
488 lambda,&dvdlambda,fr->ewald_table);
489 PRINT_SEPDVDL("Ewald long-range",Vlr,dvdlambda);
493 gmx_fatal(FARGS,"No such electrostatics method implemented %s",
494 eel_names[fr->eeltype]);
498 gmx_fatal(FARGS,"Error %d in long range electrostatics routine %s",
499 status,EELTYPE(fr->eeltype));
501 enerd->dvdl_lin += dvdlambda;
502 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
505 fprintf(debug,"Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
506 Vlr,Vcorr,enerd->term[F_COUL_RECIP]);
507 pr_rvecs(debug,0,"vir_el_recip after corr",fr->vir_el_recip,DIM);
508 pr_rvecs(debug,0,"fshift after LR Corrections",fr->fshift,SHIFTS);
513 if (EEL_RF(fr->eeltype))
517 if (fr->eeltype != eelRF_NEC)
519 enerd->term[F_RF_EXCL] =
520 RF_excl_correction(fplog,fr,graph,md,excl,x,f,
521 fr->fshift,&pbc,lambda,&dvdlambda);
524 enerd->dvdl_lin += dvdlambda;
525 PRINT_SEPDVDL("RF exclusion correction",
526 enerd->term[F_RF_EXCL],dvdlambda);
534 print_nrnb(debug,nrnb);
542 MPI_Barrier(cr->mpi_comm_mygroup);
545 if (fr->timesteps == 11)
547 fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
548 cr->nodeid, gmx_step_str(fr->timesteps,buf),
549 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
550 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
558 pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
563 void init_enerdata(int ngener,int n_flambda,gmx_enerdata_t *enerd)
567 for(i=0; i<F_NRE; i++)
575 fprintf(debug,"Creating %d sized group matrix for energies\n",n2);
577 enerd->grpp.nener = n2;
578 for(i=0; (i<egNR); i++)
580 snew(enerd->grpp.ener[i],n2);
585 enerd->n_lambda = 1 + n_flambda;
586 snew(enerd->enerpart_lambda,enerd->n_lambda);
594 void destroy_enerdata(gmx_enerdata_t *enerd)
598 for(i=0; (i<egNR); i++)
600 sfree(enerd->grpp.ener[i]);
605 sfree(enerd->enerpart_lambda);
609 static real sum_v(int n,real v[])
621 void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd)
623 gmx_grppairener_t *grpp;
630 /* Accumulate energies */
631 epot[F_COUL_SR] = sum_v(grpp->nener,grpp->ener[egCOULSR]);
632 epot[F_LJ] = sum_v(grpp->nener,grpp->ener[egLJSR]);
633 epot[F_LJ14] = sum_v(grpp->nener,grpp->ener[egLJ14]);
634 epot[F_COUL14] = sum_v(grpp->nener,grpp->ener[egCOUL14]);
635 epot[F_COUL_LR] = sum_v(grpp->nener,grpp->ener[egCOULLR]);
636 epot[F_LJ_LR] = sum_v(grpp->nener,grpp->ener[egLJLR]);
637 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
638 epot[F_GBPOL] += sum_v(grpp->nener,grpp->ener[egGB]);
640 /* lattice part of LR doesnt belong to any group
641 * and has been added earlier
643 epot[F_BHAM] = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
644 epot[F_BHAM_LR] = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
647 for(i=0; (i<F_EPOT); i++)
648 if (i != F_DISRESVIOL && i != F_ORIRESDEV && i != F_DIHRESVIOL)
649 epot[F_EPOT] += epot[i];
652 void sum_dhdl(gmx_enerdata_t *enerd,double lambda,t_inputrec *ir)
655 double dlam,dhdl_lin;
657 enerd->term[F_DVDL] = enerd->dvdl_lin + enerd->dvdl_nonlin;
661 fprintf(debug,"dvdl: %f, non-linear %f + linear %f\n",
662 enerd->term[F_DVDL],enerd->dvdl_nonlin,enerd->dvdl_lin);
665 /* Notes on the foreign lambda free energy difference evaluation:
666 * Adding the potential and ekin terms that depend linearly on lambda
667 * as delta lam * dvdl to the energy differences is exact.
668 * For the constraint dvdl this is not exact, but we have no other option.
669 * For the non-bonded LR term we assume that the soft-core (if present)
670 * no longer affects the energy beyond the short-range cut-off,
671 * which is a very good approximation (except for exotic settings).
673 for(i=1; i<enerd->n_lambda; i++)
675 dlam = (ir->flambda[i-1] - lambda);
677 enerd->dvdl_lin + enerd->term[F_DKDL] + enerd->term[F_DHDL_CON];
680 fprintf(debug,"enerdiff lam %g: non-linear %f linear %f*%f\n",
682 enerd->enerpart_lambda[i] - enerd->enerpart_lambda[0],
685 enerd->enerpart_lambda[i] += dlam*dhdl_lin;
690 void reset_enerdata(t_grpopts *opts,
691 t_forcerec *fr,gmx_bool bNS,
692 gmx_enerdata_t *enerd,
698 /* First reset all energy components, except for the long range terms
699 * on the master at non neighbor search steps, since the long range
700 * terms have already been summed at the last neighbor search step.
702 bKeepLR = (fr->bTwinRange && !bNS);
703 for(i=0; (i<egNR); i++) {
704 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
705 for(j=0; (j<enerd->grpp.nener); j++)
706 enerd->grpp.ener[i][j] = 0.0;
709 enerd->dvdl_lin = 0.0;
710 enerd->dvdl_nonlin = 0.0;
712 /* Normal potential energy components */
713 for(i=0; (i<=F_EPOT); i++) {
714 enerd->term[i] = 0.0;
716 /* Initialize the dVdlambda term with the long range contribution */
717 enerd->term[F_DVDL] = 0.0;
718 enerd->term[F_DKDL] = 0.0;
719 enerd->term[F_DHDL_CON] = 0.0;
720 if (enerd->n_lambda > 0)
722 for(i=0; i<enerd->n_lambda; i++)
724 enerd->enerpart_lambda[i] = 0.0;