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,
138 gmx_bool bDoEpot,bSepDVDL,bSB;
142 real dvdlambda,Vsr,Vlr,Vcorr=0,vdip,vcharge;
146 gmx_enerdata_t ed_lam;
151 double t0=0.0,t1,t2,t3; /* time measurement for coarse load balancing */
154 #define PRINT_SEPDVDL(s,v,dvdl) if (bSepDVDL) fprintf(fplog,sepdvdlformat,s,v,dvdl);
157 set_pbc(&pbc,fr->ePBC,box);
160 for(i=0; (i<DIM); i++)
162 box_size[i]=box[i][i];
165 bSepDVDL=(fr->bSepDVDL && do_per_step(step,ir->nstlog));
168 /* do QMMM first if requested */
171 enerd->term[F_EQM] = calculate_QMMM(cr,x,f,fr,md);
176 fprintf(fplog,"Step %s: non-bonded V and dVdl for node %d:\n",
177 gmx_step_str(step,buf),cr->nodeid);
180 /* Call the short range functions all in one go. */
185 /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
186 #define TAKETIME FALSE
189 MPI_Barrier(cr->mpi_comm_mygroup);
196 dvdlambda = do_walls(ir,fr,box,md,x,f,lambda,
197 enerd->grpp.ener[egLJSR],nrnb);
198 PRINT_SEPDVDL("Walls",0.0,dvdlambda);
199 enerd->dvdl_lin += dvdlambda;
202 /* If doing GB, reset dvda and calculate the Born radii */
203 if (ir->implicit_solvent)
205 /* wallcycle_start(wcycle,ewcGB); */
207 for(i=0;i<born->nr;i++)
214 calc_gb_rad(cr,fr,ir,top,atype,x,&(fr->gblist),born,md,nrnb);
217 /* wallcycle_stop(wcycle, ewcGB); */
222 if (flags & GMX_FORCE_FORCES)
224 donb_flags |= GMX_DONB_FORCES;
226 do_nonbonded(cr,fr,x,f,md,excl,
228 enerd->grpp.ener[egBHAMSR] :
229 enerd->grpp.ener[egLJSR],
230 enerd->grpp.ener[egCOULSR],
231 enerd->grpp.ener[egGB],box_size,nrnb,
232 lambda,&dvdlambda,-1,-1,donb_flags);
233 /* If we do foreign lambda and we have soft-core interactions
234 * we have to recalculate the (non-linear) energies contributions.
236 if (ir->n_flambda > 0 && (flags & GMX_FORCE_DHDL) && ir->sc_alpha != 0)
238 init_enerdata(mtop->groups.grps[egcENER].nr,ir->n_flambda,&ed_lam);
240 for(i=0; i<enerd->n_lambda; i++)
242 lam_i = (i==0 ? lambda : ir->flambda[i-1]);
244 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
245 do_nonbonded(cr,fr,x,f,md,excl,
247 ed_lam.grpp.ener[egBHAMSR] :
248 ed_lam.grpp.ener[egLJSR],
249 ed_lam.grpp.ener[egCOULSR],
250 enerd->grpp.ener[egGB], box_size,nrnb,
251 lam_i,&dvdl_dum,-1,-1,
252 GMX_DONB_FOREIGNLAMBDA);
253 sum_epot(&ir->opts,&ed_lam);
254 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
256 destroy_enerdata(&ed_lam);
260 /* If we are doing GB, calculate bonded forces and apply corrections
261 * to the solvation forces */
262 if (ir->implicit_solvent) {
263 calc_gb_forces(cr,md,born,top,atype,x,f,fr,idef,
264 ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&pbc,graph,enerd);
275 if (ir->sc_alpha != 0)
277 enerd->dvdl_nonlin += dvdlambda;
281 enerd->dvdl_lin += dvdlambda;
286 for(i=0; i<enerd->grpp.nener; i++)
290 enerd->grpp.ener[egBHAMSR][i] :
291 enerd->grpp.ener[egLJSR][i])
292 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
295 PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",Vsr,dvdlambda);
300 pr_rvecs(debug,0,"fshift after SR",fr->fshift,SHIFTS);
303 /* Shift the coordinates. Must be done before bonded forces and PPPM,
304 * but is also necessary for SHAKE and update, therefore it can NOT
305 * go when no bonded forces have to be evaluated.
308 /* Here sometimes we would not need to shift with NBFonly,
309 * but we do so anyhow for consistency of the returned coordinates.
313 shift_self(graph,box,x);
316 inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
320 inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
323 /* Check whether we need to do bondeds or correct for exclusions */
325 ((flags & GMX_FORCE_BONDED)
326 || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
328 /* Since all atoms are in the rectangular or triclinic unit-cell,
329 * only single box vector shifts (2 in x) are required.
331 set_pbc_dd(&pbc,fr->ePBC,cr->dd,TRUE,box);
335 if (flags & GMX_FORCE_BONDED)
337 calc_bonds(fplog,cr->ms,
338 idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
339 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
340 fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
342 /* Check if we have to determine energy differences
343 * at foreign lambda's.
345 if (ir->n_flambda > 0 && (flags & GMX_FORCE_DHDL) &&
346 idef->ilsort != ilsortNO_FE)
348 if (idef->ilsort != ilsortFE_SORTED)
350 gmx_incons("The bonded interactions are not sorted for free energy");
352 init_enerdata(mtop->groups.grps[egcENER].nr,ir->n_flambda,&ed_lam);
354 for(i=0; i<enerd->n_lambda; i++)
356 lam_i = (i==0 ? lambda : ir->flambda[i-1]);
358 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
359 calc_bonds_lambda(fplog,
360 idef,x,fr,&pbc,graph,&ed_lam,nrnb,lam_i,md,
362 DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
363 sum_epot(&ir->opts,&ed_lam);
364 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
366 destroy_enerdata(&ed_lam);
374 if (EEL_FULL(fr->eeltype))
376 bSB = (ir->nwall == 2);
380 svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
381 box_size[ZZ] *= ir->wall_ewald_zfac;
384 clear_mat(fr->vir_el_recip);
391 Vcorr = ewald_LRcorrection(fplog,md->start,md->start+md->homenr,
394 md->nChargePerturbed ? md->chargeB : NULL,
395 excl,x,bSB ? boxs : box,mu_tot,
398 lambda,&dvdlambda,&vdip,&vcharge);
399 PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdlambda);
400 enerd->dvdl_lin += dvdlambda;
404 if (ir->ewald_geometry != eewg3D || ir->epsilon_surface != 0)
406 gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
408 /* The TPI molecule does not have exclusions with the rest
409 * of the system and no intra-molecular PME grid contributions
410 * will be calculated in gmx_pme_calc_energy.
423 case eelPMEUSERSWITCH:
425 if (cr->duty & DUTY_PME)
427 assert(fr->n_tpi >= 0);
428 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
430 pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
431 if (flags & GMX_FORCE_FORCES)
433 pme_flags |= GMX_PME_CALC_F;
435 if (flags & GMX_FORCE_VIRIAL)
437 pme_flags |= GMX_PME_CALC_ENER_VIR;
441 /* We don't calculate f, but we do want the potential */
442 pme_flags |= GMX_PME_CALC_POT;
444 wallcycle_start(wcycle,ewcPMEMESH);
445 status = gmx_pme_do(fr->pmedata,
446 md->start,md->homenr - fr->n_tpi,
448 md->chargeA,md->chargeB,
450 DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
451 DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
453 fr->vir_el_recip,fr->ewaldcoeff,
454 &Vlr,lambda,&dvdlambda,
456 *cycles_pme = wallcycle_stop(wcycle,ewcPMEMESH);
458 /* We should try to do as little computation after
459 * this as possible, because parallel PME synchronizes
460 * the nodes, so we want all load imbalance of the rest
461 * of the force calculation to be before the PME call.
462 * DD load balancing is done on the whole time of
463 * the force call (without PME).
468 /* Determine the PME grid energy of the test molecule
469 * with the PME grid potential of the other charges.
471 gmx_pme_calc_energy(fr->pmedata,fr->n_tpi,
472 x + md->homenr - fr->n_tpi,
473 md->chargeA + md->homenr - fr->n_tpi,
476 PRINT_SEPDVDL("PME mesh",Vlr,dvdlambda);
480 /* Energies and virial are obtained later from the PME nodes */
481 /* but values have to be zeroed out here */
486 Vlr = do_ewald(fplog,FALSE,ir,x,fr->f_novirsum,
487 md->chargeA,md->chargeB,
488 box_size,cr,md->homenr,
489 fr->vir_el_recip,fr->ewaldcoeff,
490 lambda,&dvdlambda,fr->ewald_table);
491 PRINT_SEPDVDL("Ewald long-range",Vlr,dvdlambda);
495 gmx_fatal(FARGS,"No such electrostatics method implemented %s",
496 eel_names[fr->eeltype]);
500 gmx_fatal(FARGS,"Error %d in long range electrostatics routine %s",
501 status,EELTYPE(fr->eeltype));
503 enerd->dvdl_lin += dvdlambda;
504 enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
507 fprintf(debug,"Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
508 Vlr,Vcorr,enerd->term[F_COUL_RECIP]);
509 pr_rvecs(debug,0,"vir_el_recip after corr",fr->vir_el_recip,DIM);
510 pr_rvecs(debug,0,"fshift after LR Corrections",fr->fshift,SHIFTS);
515 if (EEL_RF(fr->eeltype))
519 if (fr->eeltype != eelRF_NEC)
521 enerd->term[F_RF_EXCL] =
522 RF_excl_correction(fplog,fr,graph,md,excl,x,f,
523 fr->fshift,&pbc,lambda,&dvdlambda);
526 enerd->dvdl_lin += dvdlambda;
527 PRINT_SEPDVDL("RF exclusion correction",
528 enerd->term[F_RF_EXCL],dvdlambda);
536 print_nrnb(debug,nrnb);
544 MPI_Barrier(cr->mpi_comm_mygroup);
547 if (fr->timesteps == 11)
549 fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
550 cr->nodeid, gmx_step_str(fr->timesteps,buf),
551 100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
552 (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
560 pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
565 void init_enerdata(int ngener,int n_flambda,gmx_enerdata_t *enerd)
569 for(i=0; i<F_NRE; i++)
577 fprintf(debug,"Creating %d sized group matrix for energies\n",n2);
579 enerd->grpp.nener = n2;
580 for(i=0; (i<egNR); i++)
582 snew(enerd->grpp.ener[i],n2);
587 enerd->n_lambda = 1 + n_flambda;
588 snew(enerd->enerpart_lambda,enerd->n_lambda);
596 void destroy_enerdata(gmx_enerdata_t *enerd)
600 for(i=0; (i<egNR); i++)
602 sfree(enerd->grpp.ener[i]);
607 sfree(enerd->enerpart_lambda);
611 static real sum_v(int n,real v[])
623 void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd)
625 gmx_grppairener_t *grpp;
632 /* Accumulate energies */
633 epot[F_COUL_SR] = sum_v(grpp->nener,grpp->ener[egCOULSR]);
634 epot[F_LJ] = sum_v(grpp->nener,grpp->ener[egLJSR]);
635 epot[F_LJ14] = sum_v(grpp->nener,grpp->ener[egLJ14]);
636 epot[F_COUL14] = sum_v(grpp->nener,grpp->ener[egCOUL14]);
637 epot[F_COUL_LR] = sum_v(grpp->nener,grpp->ener[egCOULLR]);
638 epot[F_LJ_LR] = sum_v(grpp->nener,grpp->ener[egLJLR]);
639 /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
640 epot[F_GBPOL] += sum_v(grpp->nener,grpp->ener[egGB]);
642 /* lattice part of LR doesnt belong to any group
643 * and has been added earlier
645 epot[F_BHAM] = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
646 epot[F_BHAM_LR] = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
649 for(i=0; (i<F_EPOT); i++)
650 if (i != F_DISRESVIOL && i != F_ORIRESDEV && i != F_DIHRESVIOL)
651 epot[F_EPOT] += epot[i];
654 void sum_dhdl(gmx_enerdata_t *enerd,double lambda,t_inputrec *ir)
657 double dlam,dhdl_lin;
659 enerd->term[F_DVDL] = enerd->dvdl_lin + enerd->dvdl_nonlin;
663 fprintf(debug,"dvdl: %f, non-linear %f + linear %f\n",
664 enerd->term[F_DVDL],enerd->dvdl_nonlin,enerd->dvdl_lin);
667 /* Notes on the foreign lambda free energy difference evaluation:
668 * Adding the potential and ekin terms that depend linearly on lambda
669 * as delta lam * dvdl to the energy differences is exact.
670 * For the constraint dvdl this is not exact, but we have no other option.
671 * For the non-bonded LR term we assume that the soft-core (if present)
672 * no longer affects the energy beyond the short-range cut-off,
673 * which is a very good approximation (except for exotic settings).
675 for(i=1; i<enerd->n_lambda; i++)
677 dlam = (ir->flambda[i-1] - lambda);
679 enerd->dvdl_lin + enerd->term[F_DKDL] + enerd->term[F_DHDL_CON];
682 fprintf(debug,"enerdiff lam %g: non-linear %f linear %f*%f\n",
684 enerd->enerpart_lambda[i] - enerd->enerpart_lambda[0],
687 enerd->enerpart_lambda[i] += dlam*dhdl_lin;
692 void reset_enerdata(t_grpopts *opts,
693 t_forcerec *fr,gmx_bool bNS,
694 gmx_enerdata_t *enerd,
700 /* First reset all energy components, except for the long range terms
701 * on the master at non neighbor search steps, since the long range
702 * terms have already been summed at the last neighbor search step.
704 bKeepLR = (fr->bTwinRange && !bNS);
705 for(i=0; (i<egNR); i++) {
706 if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
707 for(j=0; (j<enerd->grpp.nener); j++)
708 enerd->grpp.ener[i][j] = 0.0;
711 enerd->dvdl_lin = 0.0;
712 enerd->dvdl_nonlin = 0.0;
714 /* Normal potential energy components */
715 for(i=0; (i<=F_EPOT); i++) {
716 enerd->term[i] = 0.0;
718 /* Initialize the dVdlambda term with the long range contribution */
719 enerd->term[F_DVDL] = 0.0;
720 enerd->term[F_DKDL] = 0.0;
721 enerd->term[F_DHDL_CON] = 0.0;
722 if (enerd->n_lambda > 0)
724 for(i=0; i<enerd->n_lambda; i++)
726 enerd->enerpart_lambda[i] = 0.0;