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