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