Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / force.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 3.2.0
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.
15
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.
20  *
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.
27  *
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.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * GROwing Monsters And Cloning Shrimps
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <math.h>
41 #include <string.h>
42 #include <assert.h>
43 #include "sysstuff.h"
44 #include "typedefs.h"
45 #include "macros.h"
46 #include "smalloc.h"
47 #include "macros.h"
48 #include "physics.h"
49 #include "force.h"
50 #include "nonbonded.h"
51 #include "names.h"
52 #include "network.h"
53 #include "pbc.h"
54 #include "ns.h"
55 #include "nrnb.h"
56 #include "bondf.h"
57 #include "mshift.h"
58 #include "txtdump.h"
59 #include "coulomb.h"
60 #include "pme.h"
61 #include "mdrun.h"
62 #include "domdec.h"
63 #include "partdec.h"
64 #include "qmmm.h"
65 #include "gmx_omp_nthreads.h"
66
67
68 void ns(FILE *fp,
69         t_forcerec *fr,
70         rvec       x[],
71         matrix     box,
72         gmx_groups_t *groups,
73         t_grpopts  *opts,
74         gmx_localtop_t *top,
75         t_mdatoms  *md,
76         t_commrec  *cr,
77         t_nrnb     *nrnb,
78         real       *lambda,
79         real       *dvdlambda,
80         gmx_grppairener_t *grppener,
81         gmx_bool       bFillGrid,
82         gmx_bool       bDoLongRangeNS)
83 {
84   char   *ptr;
85   int    nsearch;
86
87
88   if (!fr->ns.nblist_initialized)
89   {
90       init_neighbor_list(fp, fr, md->homenr);
91   }
92
93   if (fr->bTwinRange)
94     fr->nlr=0;
95
96     nsearch = search_neighbours(fp,fr,x,box,top,groups,cr,nrnb,md,
97                                 lambda,dvdlambda,grppener,
98                                 bFillGrid,bDoLongRangeNS,TRUE);
99   if (debug)
100     fprintf(debug,"nsearch = %d\n",nsearch);
101
102   /* Check whether we have to do dynamic load balancing */
103   /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
104     count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
105     &(top->idef),opts->ngener);
106   */
107   if (fr->ns.dump_nl > 0)
108     dump_nblist(fp,cr,fr,fr->ns.dump_nl);
109 }
110
111 static void reduce_thread_forces(int n,rvec *f,
112                                  tensor vir,
113                                  real *Vcorr,
114                                  int efpt_ind,real *dvdl,
115                                  int nthreads,f_thread_t *f_t)
116 {
117     int t,i;
118
119     /* This reduction can run over any number of threads */
120 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
121     for(i=0; i<n; i++)
122     {
123         for(t=1; t<nthreads; t++)
124         {
125             rvec_inc(f[i],f_t[t].f[i]);
126         }
127     }
128     for(t=1; t<nthreads; t++)
129     {
130         *Vcorr += f_t[t].Vcorr;
131         *dvdl  += f_t[t].dvdl[efpt_ind];
132         m_add(vir,f_t[t].vir,vir);
133     }
134 }
135
136 void do_force_lowlevel(FILE       *fplog,   gmx_large_int_t step,
137                        t_forcerec *fr,      t_inputrec *ir,
138                        t_idef     *idef,    t_commrec  *cr,
139                        t_nrnb     *nrnb,    gmx_wallcycle_t wcycle,
140                        t_mdatoms  *md,
141                        t_grpopts  *opts,
142                        rvec       x[],      history_t  *hist,
143                        rvec       f[],
144                        rvec       f_longrange[],
145                        gmx_enerdata_t *enerd,
146                        t_fcdata   *fcd,
147                        gmx_mtop_t     *mtop,
148                        gmx_localtop_t *top,
149                        gmx_genborn_t *born,
150                        t_atomtypes *atype,
151                        gmx_bool       bBornRadii,
152                        matrix     box,
153                        t_lambda   *fepvals,
154                        real       *lambda,
155                        t_graph    *graph,
156                        t_blocka   *excl,
157                        rvec       mu_tot[],
158                        int        flags,
159                        float      *cycles_pme)
160 {
161     int     i,j,status;
162     int     donb_flags;
163     gmx_bool    bDoEpot,bSepDVDL,bSB;
164     int     pme_flags;
165     matrix  boxs;
166     rvec    box_size;
167     real    Vsr,Vlr,Vcorr=0;
168     t_pbc   pbc;
169     real    dvdgb;
170     char    buf[22];
171     gmx_enerdata_t ed_lam;
172     double  clam_i,vlam_i;
173     real    dvdl_dum[efptNR], dvdl, dvdl_nb[efptNR], lam_i[efptNR];
174     real    dvdlsum;
175
176 #ifdef GMX_MPI
177     double  t0=0.0,t1,t2,t3; /* time measurement for coarse load balancing */
178 #endif
179
180 #define PRINT_SEPDVDL(s,v,dvdlambda) if (bSepDVDL) fprintf(fplog,sepdvdlformat,s,v,dvdlambda);
181
182
183     set_pbc(&pbc,fr->ePBC,box);
184
185     /* reset free energy components */
186     for (i=0;i<efptNR;i++)
187     {
188         dvdl_nb[i]  = 0;
189         dvdl_dum[i] = 0;
190     }
191
192     /* Reset box */
193     for(i=0; (i<DIM); i++)
194     {
195         box_size[i]=box[i][i];
196     }
197
198     bSepDVDL=(fr->bSepDVDL && do_per_step(step,ir->nstlog));
199     debug_gmx();
200
201     /* do QMMM first if requested */
202     if(fr->bQMMM)
203     {
204         enerd->term[F_EQM] = calculate_QMMM(cr,x,f,fr,md);
205     }
206
207     if (bSepDVDL)
208     {
209         fprintf(fplog,"Step %s: non-bonded V and dVdl for node %d:\n",
210                 gmx_step_str(step,buf),cr->nodeid);
211     }
212
213     /* Call the short range functions all in one go. */
214
215 #ifdef GMX_MPI
216     /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
217 #define TAKETIME FALSE
218     if (TAKETIME)
219     {
220         MPI_Barrier(cr->mpi_comm_mygroup);
221         t0=MPI_Wtime();
222     }
223 #endif
224
225     if (ir->nwall)
226     {
227         /* foreign lambda component for walls */
228         dvdl = do_walls(ir,fr,box,md,x,f,lambda[efptVDW],
229                         enerd->grpp.ener[egLJSR],nrnb);
230         PRINT_SEPDVDL("Walls",0.0,dvdl);
231         enerd->dvdl_lin[efptVDW] += dvdl;
232     }
233
234         /* If doing GB, reset dvda and calculate the Born radii */
235         if (ir->implicit_solvent)
236         {
237         wallcycle_sub_start(wcycle, ewcsNONBONDED);
238
239                 for(i=0;i<born->nr;i++)
240                 {
241                         fr->dvda[i]=0;
242                 }
243
244                 if(bBornRadii)
245                 {
246                         calc_gb_rad(cr,fr,ir,top,atype,x,&(fr->gblist),born,md,nrnb);
247                 }
248
249         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
250         }
251
252     where();
253     /* We only do non-bonded calculation with group scheme here, the verlet
254      * calls are done from do_force_cutsVERLET(). */
255     if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
256     {
257         donb_flags = 0;
258         /* Add short-range interactions */
259         donb_flags |= GMX_NONBONDED_DO_SR;
260
261         if (flags & GMX_FORCE_FORCES)
262         {
263             donb_flags |= GMX_NONBONDED_DO_FORCE;
264         }
265         if (flags & GMX_FORCE_ENERGY)
266         {
267             donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
268         }
269         if (flags & GMX_FORCE_DO_LR)
270         {
271             donb_flags |= GMX_NONBONDED_DO_LR;
272         }
273
274         wallcycle_sub_start(wcycle, ewcsNONBONDED);
275         do_nonbonded(cr,fr,x,f,f_longrange,md,excl,
276                     &enerd->grpp,box_size,nrnb,
277                     lambda,dvdl_nb,-1,-1,donb_flags);
278         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
279     }
280
281     /* If we do foreign lambda and we have soft-core interactions
282      * we have to recalculate the (non-linear) energies contributions.
283      */
284     if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
285     {
286         wallcycle_sub_start(wcycle, ewcsNONBONDED);
287         init_enerdata(mtop->groups.grps[egcENER].nr,fepvals->n_lambda,&ed_lam);
288
289         for(i=0; i<enerd->n_lambda; i++)
290         {
291             for (j=0;j<efptNR;j++)
292             {
293                 lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
294             }
295             reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
296             do_nonbonded(cr,fr,x,f,f_longrange,md,excl,
297                          &(ed_lam.grpp), box_size,nrnb,
298                          lam_i,dvdl_dum,-1,-1,
299                          GMX_NONBONDED_DO_FOREIGNLAMBDA | GMX_NONBONDED_DO_SR);
300             sum_epot(&ir->opts,&ed_lam);
301             enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
302         }
303         destroy_enerdata(&ed_lam);
304         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
305     }
306     where();
307
308         /* If we are doing GB, calculate bonded forces and apply corrections
309          * to the solvation forces */
310     /* MRS: Eventually, many need to include free energy contribution here! */
311         if (ir->implicit_solvent)
312     {
313         wallcycle_sub_start(wcycle, ewcsBONDED);
314                 calc_gb_forces(cr,md,born,top,atype,x,f,fr,idef,
315                        ir->gb_algorithm,ir->sa_algorithm,nrnb,bBornRadii,&pbc,graph,enerd);
316         wallcycle_sub_stop(wcycle, ewcsBONDED);
317     }
318
319 #ifdef GMX_MPI
320     if (TAKETIME)
321     {
322         t1=MPI_Wtime();
323         fr->t_fnbf += t1-t0;
324     }
325 #endif
326
327     if (fepvals->sc_alpha!=0)
328     {
329         enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
330     }
331     else
332     {
333         enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
334     }
335
336     if (fepvals->sc_alpha!=0)
337
338         /* even though coulomb part is linear, we already added it, beacuse we
339            need to go through the vdw calculation anyway */
340     {
341         enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
342     }
343     else
344     {
345         enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
346     }
347
348     Vsr = 0;
349     if (bSepDVDL)
350     {
351         for(i=0; i<enerd->grpp.nener; i++)
352         {
353             Vsr +=
354                 (fr->bBHAM ?
355                  enerd->grpp.ener[egBHAMSR][i] :
356                  enerd->grpp.ener[egLJSR][i])
357                 + enerd->grpp.ener[egCOULSR][i] + enerd->grpp.ener[egGB][i];
358         }
359         dvdlsum = dvdl_nb[efptVDW] + dvdl_nb[efptCOUL];
360         PRINT_SEPDVDL("VdW and Coulomb SR particle-p.",Vsr,dvdlsum);
361     }
362     debug_gmx();
363
364
365     if (debug)
366     {
367         pr_rvecs(debug,0,"fshift after SR",fr->fshift,SHIFTS);
368     }
369
370     /* Shift the coordinates. Must be done before bonded forces and PPPM,
371      * but is also necessary for SHAKE and update, therefore it can NOT
372      * go when no bonded forces have to be evaluated.
373      */
374
375     /* Here sometimes we would not need to shift with NBFonly,
376      * but we do so anyhow for consistency of the returned coordinates.
377      */
378     if (graph)
379     {
380         shift_self(graph,box,x);
381         if (TRICLINIC(box))
382         {
383             inc_nrnb(nrnb,eNR_SHIFTX,2*graph->nnodes);
384         }
385         else
386         {
387             inc_nrnb(nrnb,eNR_SHIFTX,graph->nnodes);
388         }
389     }
390     /* Check whether we need to do bondeds or correct for exclusions */
391     if (fr->bMolPBC &&
392         ((flags & GMX_FORCE_BONDED)
393          || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype)))
394     {
395         /* Since all atoms are in the rectangular or triclinic unit-cell,
396          * only single box vector shifts (2 in x) are required.
397          */
398         set_pbc_dd(&pbc,fr->ePBC,cr->dd,TRUE,box);
399     }
400     debug_gmx();
401
402     if (flags & GMX_FORCE_BONDED)
403     {
404         wallcycle_sub_start(wcycle, ewcsBONDED);
405         calc_bonds(fplog,cr->ms,
406                    idef,x,hist,f,fr,&pbc,graph,enerd,nrnb,lambda,md,fcd,
407                    DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL, atype, born,
408                    flags,
409                    fr->bSepDVDL && do_per_step(step,ir->nstlog),step);
410
411         /* Check if we have to determine energy differences
412          * at foreign lambda's.
413          */
414         if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
415             idef->ilsort != ilsortNO_FE)
416         {
417             if (idef->ilsort != ilsortFE_SORTED)
418             {
419                 gmx_incons("The bonded interactions are not sorted for free energy");
420             }
421             init_enerdata(mtop->groups.grps[egcENER].nr,fepvals->n_lambda,&ed_lam);
422
423             for(i=0; i<enerd->n_lambda; i++)
424             {
425                 reset_enerdata(&ir->opts,fr,TRUE,&ed_lam,FALSE);
426                 for (j=0;j<efptNR;j++)
427                 {
428                     lam_i[j] = (i==0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
429                 }
430                 calc_bonds_lambda(fplog,idef,x,fr,&pbc,graph,&ed_lam,nrnb,lam_i,md,
431                                   fcd,DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
432                 sum_epot(&ir->opts,&ed_lam);
433                 enerd->enerpart_lambda[i] += ed_lam.term[F_EPOT];
434             }
435             destroy_enerdata(&ed_lam);
436         }
437         debug_gmx();
438
439         wallcycle_sub_stop(wcycle, ewcsBONDED);
440     }
441
442     where();
443
444     *cycles_pme = 0;
445     if (EEL_FULL(fr->eeltype))
446     {
447         bSB = (ir->nwall == 2);
448         if (bSB)
449         {
450             copy_mat(box,boxs);
451             svmul(ir->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
452             box_size[ZZ] *= ir->wall_ewald_zfac;
453         }
454
455         clear_mat(fr->vir_el_recip);
456
457         if (fr->bEwald)
458         {
459             Vcorr = 0;
460             dvdl  = 0;
461
462             /* With the Verlet scheme exclusion forces are calculated
463              * in the non-bonded kernel.
464              */
465             /* The TPI molecule does not have exclusions with the rest
466              * of the system and no intra-molecular PME grid contributions
467              * will be calculated in gmx_pme_calc_energy.
468              */
469             if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
470                 ir->ewald_geometry != eewg3D ||
471                 ir->epsilon_surface != 0)
472             {
473                 int nthreads,t;
474
475                 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
476
477                 if (fr->n_tpi > 0)
478                 {
479                     gmx_fatal(FARGS,"TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
480                 }
481
482                 nthreads = gmx_omp_nthreads_get(emntBonded);
483 #pragma omp parallel for num_threads(nthreads) schedule(static)
484                 for(t=0; t<nthreads; t++)
485                 {
486                     int s,e,i;
487                     rvec *fnv;
488                     tensor *vir;
489                     real *Vcorrt,*dvdlt;
490                     if (t == 0)
491                     {
492                         fnv    = fr->f_novirsum;
493                         vir    = &fr->vir_el_recip;
494                         Vcorrt = &Vcorr;
495                         dvdlt  = &dvdl;
496                     }
497                     else
498                     {
499                         fnv    = fr->f_t[t].f;
500                         vir    = &fr->f_t[t].vir;
501                         Vcorrt = &fr->f_t[t].Vcorr;
502                         dvdlt  = &fr->f_t[t].dvdl[efptCOUL];
503                         for(i=0; i<fr->natoms_force; i++)
504                         {
505                             clear_rvec(fnv[i]);
506                         }
507                         clear_mat(*vir);
508                     }
509                     *dvdlt = 0;
510                     *Vcorrt =
511                         ewald_LRcorrection(fplog,
512                                            fr->excl_load[t],fr->excl_load[t+1],
513                                            cr,t,fr,
514                                            md->chargeA,
515                                            md->nChargePerturbed ? md->chargeB : NULL,
516                                            ir->cutoff_scheme != ecutsVERLET,
517                                            excl,x,bSB ? boxs : box,mu_tot,
518                                            ir->ewald_geometry,
519                                            ir->epsilon_surface,
520                                            fnv,*vir,
521                                            lambda[efptCOUL],dvdlt);
522                 }
523                 if (nthreads > 1)
524                 {
525                     reduce_thread_forces(fr->natoms_force,fr->f_novirsum,
526                                          fr->vir_el_recip,
527                                          &Vcorr,efptCOUL,&dvdl,
528                                          nthreads,fr->f_t);
529                 }
530
531                 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
532             }
533
534             if (fr->n_tpi == 0)
535             {
536                 Vcorr += ewald_charge_correction(cr,fr,lambda[efptCOUL],box,
537                                                  &dvdl,fr->vir_el_recip);
538             }
539
540             PRINT_SEPDVDL("Ewald excl./charge/dip. corr.",Vcorr,dvdl);
541             enerd->dvdl_lin[efptCOUL] += dvdl;
542         }
543
544         status = 0;
545         Vlr  = 0;
546         dvdl = 0;
547         switch (fr->eeltype)
548         {
549         case eelPME:
550         case eelPMESWITCH:
551         case eelPMEUSER:
552         case eelPMEUSERSWITCH:
553         case eelP3M_AD:
554             if (cr->duty & DUTY_PME)
555             {
556                 assert(fr->n_tpi >= 0);
557                 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
558                 {
559                     pme_flags = GMX_PME_SPREAD_Q | GMX_PME_SOLVE;
560                     if (flags & GMX_FORCE_FORCES)
561                     {
562                         pme_flags |= GMX_PME_CALC_F;
563                     }
564                     if (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY))
565                     {
566                         pme_flags |= GMX_PME_CALC_ENER_VIR;
567                     }
568                     if (fr->n_tpi > 0)
569                     {
570                         /* We don't calculate f, but we do want the potential */
571                         pme_flags |= GMX_PME_CALC_POT;
572                     }
573                     wallcycle_start(wcycle,ewcPMEMESH);
574                     status = gmx_pme_do(fr->pmedata,
575                                         md->start,md->homenr - fr->n_tpi,
576                                         x,fr->f_novirsum,
577                                         md->chargeA,md->chargeB,
578                                         bSB ? boxs : box,cr,
579                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
580                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
581                                         nrnb,wcycle,
582                                         fr->vir_el_recip,fr->ewaldcoeff,
583                                         &Vlr,lambda[efptCOUL],&dvdl,
584                                         pme_flags);
585                     *cycles_pme = wallcycle_stop(wcycle,ewcPMEMESH);
586
587                     /* We should try to do as little computation after
588                      * this as possible, because parallel PME synchronizes
589                      * the nodes, so we want all load imbalance of the rest
590                      * of the force calculation to be before the PME call.
591                      * DD load balancing is done on the whole time of
592                      * the force call (without PME).
593                      */
594                 }
595                 if (fr->n_tpi > 0)
596                 {
597                     /* Determine the PME grid energy of the test molecule
598                      * with the PME grid potential of the other charges.
599                      */
600                     gmx_pme_calc_energy(fr->pmedata,fr->n_tpi,
601                                         x + md->homenr - fr->n_tpi,
602                                         md->chargeA + md->homenr - fr->n_tpi,
603                                         &Vlr);
604                 }
605                 PRINT_SEPDVDL("PME mesh",Vlr,dvdl);
606             }
607             break;
608         case eelEWALD:
609             Vlr = do_ewald(fplog,FALSE,ir,x,fr->f_novirsum,
610                            md->chargeA,md->chargeB,
611                            box_size,cr,md->homenr,
612                            fr->vir_el_recip,fr->ewaldcoeff,
613                            lambda[efptCOUL],&dvdl,fr->ewald_table);
614             PRINT_SEPDVDL("Ewald long-range",Vlr,dvdl);
615             break;
616         default:
617             gmx_fatal(FARGS,"No such electrostatics method implemented %s",
618                       eel_names[fr->eeltype]);
619         }
620         if (status != 0)
621         {
622             gmx_fatal(FARGS,"Error %d in long range electrostatics routine %s",
623                       status,EELTYPE(fr->eeltype));
624                 }
625         /* Note that with separate PME nodes we get the real energies later */
626         enerd->dvdl_lin[efptCOUL] += dvdl;
627         enerd->term[F_COUL_RECIP] = Vlr + Vcorr;
628         if (debug)
629         {
630             fprintf(debug,"Vlr = %g, Vcorr = %g, Vlr_corr = %g\n",
631                     Vlr,Vcorr,enerd->term[F_COUL_RECIP]);
632             pr_rvecs(debug,0,"vir_el_recip after corr",fr->vir_el_recip,DIM);
633             pr_rvecs(debug,0,"fshift after LR Corrections",fr->fshift,SHIFTS);
634         }
635     }
636     else
637     {
638         if (EEL_RF(fr->eeltype))
639         {
640             /* With the Verlet scheme exclusion forces are calculated
641              * in the non-bonded kernel.
642              */
643             if (ir->cutoff_scheme != ecutsVERLET && fr->eeltype != eelRF_NEC)
644             {
645                 dvdl = 0;
646                 enerd->term[F_RF_EXCL] =
647                     RF_excl_correction(fplog,fr,graph,md,excl,x,f,
648                                        fr->fshift,&pbc,lambda[efptCOUL],&dvdl);
649             }
650
651             enerd->dvdl_lin[efptCOUL] += dvdl;
652             PRINT_SEPDVDL("RF exclusion correction",
653                           enerd->term[F_RF_EXCL],dvdl);
654         }
655     }
656     where();
657     debug_gmx();
658
659     if (debug)
660     {
661         print_nrnb(debug,nrnb);
662     }
663     debug_gmx();
664
665 #ifdef GMX_MPI
666     if (TAKETIME)
667     {
668         t2=MPI_Wtime();
669         MPI_Barrier(cr->mpi_comm_mygroup);
670         t3=MPI_Wtime();
671         fr->t_wait += t3-t2;
672         if (fr->timesteps == 11)
673         {
674             fprintf(stderr,"* PP load balancing info: node %d, step %s, rel wait time=%3.0f%% , load string value: %7.2f\n",
675                     cr->nodeid, gmx_step_str(fr->timesteps,buf),
676                     100*fr->t_wait/(fr->t_wait+fr->t_fnbf),
677                     (fr->t_fnbf+fr->t_wait)/fr->t_fnbf);
678         }
679         fr->timesteps++;
680     }
681 #endif
682
683     if (debug)
684     {
685         pr_rvecs(debug,0,"fshift after bondeds",fr->fshift,SHIFTS);
686     }
687
688 }
689
690 void init_enerdata(int ngener,int n_lambda,gmx_enerdata_t *enerd)
691 {
692     int i,n2;
693
694     for(i=0; i<F_NRE; i++)
695     {
696         enerd->term[i] = 0;
697     }
698
699
700     for(i=0; i<efptNR; i++) {
701         enerd->dvdl_lin[i]  = 0;
702         enerd->dvdl_nonlin[i]  = 0;
703     }
704
705     n2=ngener*ngener;
706     if (debug)
707     {
708         fprintf(debug,"Creating %d sized group matrix for energies\n",n2);
709     }
710     enerd->grpp.nener = n2;
711     for(i=0; (i<egNR); i++)
712     {
713         snew(enerd->grpp.ener[i],n2);
714     }
715
716     if (n_lambda)
717     {
718         enerd->n_lambda = 1 + n_lambda;
719         snew(enerd->enerpart_lambda,enerd->n_lambda);
720     }
721     else
722     {
723         enerd->n_lambda = 0;
724     }
725 }
726
727 void destroy_enerdata(gmx_enerdata_t *enerd)
728 {
729     int i;
730
731     for(i=0; (i<egNR); i++)
732     {
733         sfree(enerd->grpp.ener[i]);
734     }
735
736     if (enerd->n_lambda)
737     {
738         sfree(enerd->enerpart_lambda);
739     }
740 }
741
742 static real sum_v(int n,real v[])
743 {
744   real t;
745   int  i;
746
747   t = 0.0;
748   for(i=0; (i<n); i++)
749     t = t + v[i];
750
751   return t;
752 }
753
754 void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd)
755 {
756   gmx_grppairener_t *grpp;
757   real *epot;
758   int i;
759
760   grpp = &enerd->grpp;
761   epot = enerd->term;
762
763   /* Accumulate energies */
764   epot[F_COUL_SR]  = sum_v(grpp->nener,grpp->ener[egCOULSR]);
765   epot[F_LJ]       = sum_v(grpp->nener,grpp->ener[egLJSR]);
766   epot[F_LJ14]     = sum_v(grpp->nener,grpp->ener[egLJ14]);
767   epot[F_COUL14]   = sum_v(grpp->nener,grpp->ener[egCOUL14]);
768   epot[F_COUL_LR]  = sum_v(grpp->nener,grpp->ener[egCOULLR]);
769   epot[F_LJ_LR]    = sum_v(grpp->nener,grpp->ener[egLJLR]);
770   /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
771   epot[F_GBPOL]   += sum_v(grpp->nener,grpp->ener[egGB]);
772
773 /* lattice part of LR doesnt belong to any group
774  * and has been added earlier
775  */
776   epot[F_BHAM]     = sum_v(grpp->nener,grpp->ener[egBHAMSR]);
777   epot[F_BHAM_LR]  = sum_v(grpp->nener,grpp->ener[egBHAMLR]);
778
779   epot[F_EPOT] = 0;
780   for(i=0; (i<F_EPOT); i++)
781   {
782       if (i != F_DISRESVIOL && i != F_ORIRESDEV)
783       {
784           epot[F_EPOT] += epot[i];
785       }
786   }
787 }
788
789 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
790 {
791     int i,j,index;
792     double dlam;
793
794     enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW];  /* include dispersion correction */
795     enerd->term[F_DVDL] = 0.0;
796     for (i=0;i<efptNR;i++)
797     {
798         if (fepvals->separate_dvdl[i])
799         {
800             /* could this be done more readably/compactly? */
801             switch (i) {
802             case (efptCOUL):
803                 index = F_DVDL_COUL;
804                 break;
805             case (efptVDW):
806                 index = F_DVDL_VDW;
807                 break;
808             case (efptBONDED):
809                 index = F_DVDL_BONDED;
810                 break;
811             case (efptRESTRAINT):
812                 index = F_DVDL_RESTRAINT;
813                 break;
814             case (efptMASS):
815                 index = F_DKDL;
816                 break;
817             default:
818                 index = F_DVDL;
819                 break;
820             }
821             enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
822             if (debug)
823             {
824                 fprintf(debug,"dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
825                         efpt_names[i],i,enerd->term[index],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
826             }
827         }
828         else
829         {
830             enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
831             if (debug)
832             {
833                 fprintf(debug,"dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
834                         efpt_names[0],i,enerd->term[F_DVDL],enerd->dvdl_nonlin[i],enerd->dvdl_lin[i]);
835             }
836         }
837     }
838
839     /* Notes on the foreign lambda free energy difference evaluation:
840      * Adding the potential and ekin terms that depend linearly on lambda
841      * as delta lam * dvdl to the energy differences is exact.
842      * For the constraints this is not exact, but we have no other option
843      * without literally changing the lengths and reevaluating the energies at each step.
844      * (try to remedy this post 4.6 - MRS)
845      * For the non-bonded LR term we assume that the soft-core (if present)
846      * no longer affects the energy beyond the short-range cut-off,
847      * which is a very good approximation (except for exotic settings).
848      * (investigate how to overcome this post 4.6 - MRS)
849      */
850
851     for(i=0; i<fepvals->n_lambda; i++)
852     {                                         /* note we are iterating over fepvals here!
853                                                  For the current lam, dlam = 0 automatically,
854                                                  so we don't need to add anything to the
855                                                  enerd->enerpart_lambda[0] */
856
857         /* we don't need to worry about dvdl contributions to the current lambda, because
858            it's automatically zero */
859
860         /* first kinetic energy term */
861         dlam = (fepvals->all_lambda[efptMASS][i] - lambda[efptMASS]);
862
863         enerd->enerpart_lambda[i+1] += enerd->term[F_DKDL]*dlam;
864
865         for (j=0;j<efptNR;j++)
866         {
867             if (j==efptMASS) {continue;} /* no other mass term to worry about */
868
869             dlam = (fepvals->all_lambda[j][i]-lambda[j]);
870             enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
871             if (debug)
872             {
873                 fprintf(debug,"enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
874                         fepvals->all_lambda[j][i],efpt_names[j],
875                         (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
876                         dlam,enerd->dvdl_lin[j]);
877             }
878         }
879     }
880 }
881
882 void reset_enerdata(t_grpopts *opts,
883                     t_forcerec *fr,gmx_bool bNS,
884                     gmx_enerdata_t *enerd,
885                     gmx_bool bMaster)
886 {
887     gmx_bool bKeepLR;
888     int  i,j;
889
890     /* First reset all energy components, except for the long range terms
891      * on the master at non neighbor search steps, since the long range
892      * terms have already been summed at the last neighbor search step.
893      */
894     bKeepLR = (fr->bTwinRange && !bNS);
895     for(i=0; (i<egNR); i++) {
896         if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR))) {
897             for(j=0; (j<enerd->grpp.nener); j++)
898                 enerd->grpp.ener[i][j] = 0.0;
899         }
900     }
901     for (i=0;i<efptNR;i++)
902     {
903         enerd->dvdl_lin[i]    = 0.0;
904         enerd->dvdl_nonlin[i] = 0.0;
905     }
906
907     /* Normal potential energy components */
908     for(i=0; (i<=F_EPOT); i++) {
909         enerd->term[i] = 0.0;
910     }
911     /* Initialize the dVdlambda term with the long range contribution */
912     /* Initialize the dvdl term with the long range contribution */
913     enerd->term[F_DVDL]            = 0.0;
914     enerd->term[F_DVDL_COUL]       = 0.0;
915     enerd->term[F_DVDL_VDW]        = 0.0;
916     enerd->term[F_DVDL_BONDED]     = 0.0;
917     enerd->term[F_DVDL_RESTRAINT]  = 0.0;
918     enerd->term[F_DKDL]            = 0.0;
919     if (enerd->n_lambda > 0)
920     {
921         for(i=0; i<enerd->n_lambda; i++)
922         {
923             enerd->enerpart_lambda[i] = 0.0;
924         }
925     }
926 }