Renamed bonded module as 'listed-forces'
[alexxy/gromacs.git] / src / gromacs / mdlib / force.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "gromacs/legacyheaders/force.h"
40
41 #include "config.h"
42
43 #include <assert.h>
44 #include <math.h>
45 #include <string.h>
46
47 #include "gromacs/legacyheaders/coulomb.h"
48 #include "gromacs/legacyheaders/domdec.h"
49 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/mdrun.h"
52 #include "gromacs/legacyheaders/names.h"
53 #include "gromacs/legacyheaders/network.h"
54 #include "gromacs/legacyheaders/nonbonded.h"
55 #include "gromacs/legacyheaders/nrnb.h"
56 #include "gromacs/legacyheaders/ns.h"
57 #include "gromacs/legacyheaders/pme.h"
58 #include "gromacs/legacyheaders/qmmm.h"
59 #include "gromacs/legacyheaders/txtdump.h"
60 #include "gromacs/legacyheaders/typedefs.h"
61 #include "gromacs/legacyheaders/types/commrec.h"
62 #include "gromacs/listed-forces/bonded.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/pbcutil/ishift.h"
65 #include "gromacs/pbcutil/mshift.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/smalloc.h"
70
71 void ns(FILE              *fp,
72         t_forcerec        *fr,
73         matrix             box,
74         gmx_groups_t      *groups,
75         gmx_localtop_t    *top,
76         t_mdatoms         *md,
77         t_commrec         *cr,
78         t_nrnb            *nrnb,
79         gmx_bool           bFillGrid,
80         gmx_bool           bDoLongRangeNS)
81 {
82     char   *ptr;
83     int     nsearch;
84
85
86     if (!fr->ns.nblist_initialized)
87     {
88         init_neighbor_list(fp, fr, md->homenr);
89     }
90
91     if (fr->bTwinRange)
92     {
93         fr->nlr = 0;
94     }
95
96     nsearch = search_neighbours(fp, fr, box, top, groups, cr, nrnb, md,
97                                 bFillGrid, bDoLongRangeNS);
98     if (debug)
99     {
100         fprintf(debug, "nsearch = %d\n", nsearch);
101     }
102
103     /* Check whether we have to do dynamic load balancing */
104     /*if ((nsb->nstDlb > 0) && (mod(step,nsb->nstDlb) == 0))
105        count_nb(cr,nsb,&(top->blocks[ebCGS]),nns,fr->nlr,
106        &(top->idef),opts->ngener);
107      */
108     if (fr->ns.dump_nl > 0)
109     {
110         dump_nblist(fp, cr, fr, fr->ns.dump_nl);
111     }
112 }
113
114 static void reduce_thread_forces(int n, rvec *f,
115                                  tensor vir_q, tensor vir_lj,
116                                  real *Vcorr_q, real *Vcorr_lj,
117                                  real *dvdl_q, real *dvdl_lj,
118                                  int nthreads, f_thread_t *f_t)
119 {
120     int t, i;
121
122     /* This reduction can run over any number of threads */
123 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntBonded)) private(t) schedule(static)
124     for (i = 0; i < n; i++)
125     {
126         for (t = 1; t < nthreads; t++)
127         {
128             rvec_inc(f[i], f_t[t].f[i]);
129         }
130     }
131     for (t = 1; t < nthreads; t++)
132     {
133         *Vcorr_q  += f_t[t].Vcorr_q;
134         *Vcorr_lj += f_t[t].Vcorr_lj;
135         *dvdl_q   += f_t[t].dvdl[efptCOUL];
136         *dvdl_lj  += f_t[t].dvdl[efptVDW];
137         m_add(vir_q, f_t[t].vir_q, vir_q);
138         m_add(vir_lj, f_t[t].vir_lj, vir_lj);
139     }
140 }
141
142 void do_force_lowlevel(t_forcerec *fr,      t_inputrec *ir,
143                        t_idef     *idef,    t_commrec  *cr,
144                        t_nrnb     *nrnb,    gmx_wallcycle_t wcycle,
145                        t_mdatoms  *md,
146                        rvec       x[],      history_t  *hist,
147                        rvec       f[],
148                        rvec       f_longrange[],
149                        gmx_enerdata_t *enerd,
150                        t_fcdata   *fcd,
151                        gmx_localtop_t *top,
152                        gmx_genborn_t *born,
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;
164     int         donb_flags;
165     gmx_bool    bDoEpot, bSB;
166     int         pme_flags;
167     matrix      boxs;
168     rvec        box_size;
169     t_pbc       pbc;
170     char        buf[22];
171     double      clam_i, vlam_i;
172     real        dvdl_dum[efptNR], dvdl_nb[efptNR], lam_i[efptNR];
173     real        dvdl_q, dvdl_lj;
174
175 #ifdef GMX_MPI
176     double  t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
177 #endif
178
179     set_pbc(&pbc, fr->ePBC, box);
180
181     /* reset free energy components */
182     for (i = 0; i < efptNR; i++)
183     {
184         dvdl_nb[i]  = 0;
185         dvdl_dum[i] = 0;
186     }
187
188     /* Reset box */
189     for (i = 0; (i < DIM); i++)
190     {
191         box_size[i] = box[i][i];
192     }
193
194     debug_gmx();
195
196     /* do QMMM first if requested */
197     if (fr->bQMMM)
198     {
199         enerd->term[F_EQM] = calculate_QMMM(cr, x, f, fr);
200     }
201
202     /* Call the short range functions all in one go. */
203
204 #ifdef GMX_MPI
205     /*#define TAKETIME ((cr->npmenodes) && (fr->timesteps < 12))*/
206 #define TAKETIME FALSE
207     if (TAKETIME)
208     {
209         MPI_Barrier(cr->mpi_comm_mygroup);
210         t0 = MPI_Wtime();
211     }
212 #endif
213
214     if (ir->nwall)
215     {
216         /* foreign lambda component for walls */
217         real dvdl_walls = do_walls(ir, fr, box, md, x, f, lambda[efptVDW],
218                                    enerd->grpp.ener[egLJSR], nrnb);
219         enerd->dvdl_lin[efptVDW] += dvdl_walls;
220     }
221
222     /* If doing GB, reset dvda and calculate the Born radii */
223     if (ir->implicit_solvent)
224     {
225         wallcycle_sub_start(wcycle, ewcsNONBONDED);
226
227         for (i = 0; i < born->nr; i++)
228         {
229             fr->dvda[i] = 0;
230         }
231
232         if (bBornRadii)
233         {
234             calc_gb_rad(cr, fr, ir, top, x, &(fr->gblist), born, md, nrnb);
235         }
236
237         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
238     }
239
240     where();
241     /* We only do non-bonded calculation with group scheme here, the verlet
242      * calls are done from do_force_cutsVERLET(). */
243     if (fr->cutoff_scheme == ecutsGROUP && (flags & GMX_FORCE_NONBONDED))
244     {
245         donb_flags = 0;
246         /* Add short-range interactions */
247         donb_flags |= GMX_NONBONDED_DO_SR;
248
249         /* Currently all group scheme kernels always calculate (shift-)forces */
250         if (flags & GMX_FORCE_FORCES)
251         {
252             donb_flags |= GMX_NONBONDED_DO_FORCE;
253         }
254         if (flags & GMX_FORCE_VIRIAL)
255         {
256             donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
257         }
258         if (flags & GMX_FORCE_ENERGY)
259         {
260             donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
261         }
262         if (flags & GMX_FORCE_DO_LR)
263         {
264             donb_flags |= GMX_NONBONDED_DO_LR;
265         }
266
267         wallcycle_sub_start(wcycle, ewcsNONBONDED);
268         do_nonbonded(fr, x, f, f_longrange, md, excl,
269                      &enerd->grpp, nrnb,
270                      lambda, dvdl_nb, -1, -1, donb_flags);
271
272         /* If we do foreign lambda and we have soft-core interactions
273          * we have to recalculate the (non-linear) energies contributions.
274          */
275         if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
276         {
277             for (i = 0; i < enerd->n_lambda; i++)
278             {
279                 for (j = 0; j < efptNR; j++)
280                 {
281                     lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
282                 }
283                 reset_foreign_enerdata(enerd);
284                 do_nonbonded(fr, x, f, f_longrange, md, excl,
285                              &(enerd->foreign_grpp), nrnb,
286                              lam_i, dvdl_dum, -1, -1,
287                              (donb_flags & ~GMX_NONBONDED_DO_FORCE) | GMX_NONBONDED_DO_FOREIGNLAMBDA);
288                 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
289                 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
290             }
291         }
292         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
293         where();
294     }
295
296     /* If we are doing GB, calculate bonded forces and apply corrections
297      * to the solvation forces */
298     /* MRS: Eventually, many need to include free energy contribution here! */
299     if (ir->implicit_solvent)
300     {
301         wallcycle_sub_start(wcycle, ewcsLISTED);
302         calc_gb_forces(cr, md, born, top, x, f, fr, idef,
303                        ir->gb_algorithm, ir->sa_algorithm, nrnb, &pbc, graph, enerd);
304         wallcycle_sub_stop(wcycle, ewcsLISTED);
305     }
306
307 #ifdef GMX_MPI
308     if (TAKETIME)
309     {
310         t1          = MPI_Wtime();
311         fr->t_fnbf += t1-t0;
312     }
313 #endif
314
315     if (fepvals->sc_alpha != 0)
316     {
317         enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
318     }
319     else
320     {
321         enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
322     }
323
324     if (fepvals->sc_alpha != 0)
325
326     /* even though coulomb part is linear, we already added it, beacuse we
327        need to go through the vdw calculation anyway */
328     {
329         enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
330     }
331     else
332     {
333         enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
334     }
335
336     debug_gmx();
337
338
339     if (debug)
340     {
341         pr_rvecs(debug, 0, "fshift after SR", fr->fshift, SHIFTS);
342     }
343
344     /* Shift the coordinates. Must be done before listed forces and PPPM,
345      * but is also necessary for SHAKE and update, therefore it can NOT
346      * go when no listed forces have to be evaluated.
347      */
348
349     /* Here sometimes we would not need to shift with NBFonly,
350      * but we do so anyhow for consistency of the returned coordinates.
351      */
352     if (graph)
353     {
354         shift_self(graph, box, x);
355         if (TRICLINIC(box))
356         {
357             inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
358         }
359         else
360         {
361             inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
362         }
363     }
364     /* Check whether we need to do listed interactions or correct for exclusions */
365     if (fr->bMolPBC &&
366         ((flags & GMX_FORCE_LISTED)
367          || EEL_RF(fr->eeltype) || EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype)))
368     {
369         /* Since all atoms are in the rectangular or triclinic unit-cell,
370          * only single box vector shifts (2 in x) are required.
371          */
372         set_pbc_dd(&pbc, fr->ePBC, cr->dd, TRUE, box);
373     }
374     debug_gmx();
375
376     if (flags & GMX_FORCE_LISTED)
377     {
378         wallcycle_sub_start(wcycle, ewcsLISTED);
379         calc_bonds(cr->ms,
380                    idef, (const rvec *) x, hist, f, fr, &pbc, graph, enerd, nrnb, lambda, md, fcd,
381                    DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL,
382                    flags);
383
384         /* Check if we have to determine energy differences
385          * at foreign lambda's.
386          */
387         if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) &&
388             idef->ilsort != ilsortNO_FE)
389         {
390             if (idef->ilsort != ilsortFE_SORTED)
391             {
392                 gmx_incons("The bonded interactions are not sorted for free energy");
393             }
394             for (i = 0; i < enerd->n_lambda; i++)
395             {
396                 reset_foreign_enerdata(enerd);
397                 for (j = 0; j < efptNR; j++)
398                 {
399                     lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
400                 }
401                 calc_bonds_lambda(idef, (const rvec *) x, fr, &pbc, graph, &(enerd->foreign_grpp), enerd->foreign_term, nrnb, lam_i, md,
402                                   fcd, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
403                 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
404                 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
405             }
406         }
407         debug_gmx();
408
409         wallcycle_sub_stop(wcycle, ewcsLISTED);
410     }
411
412     where();
413
414     *cycles_pme = 0;
415     clear_mat(fr->vir_el_recip);
416     clear_mat(fr->vir_lj_recip);
417
418     /* Do long-range electrostatics and/or LJ-PME, including related short-range
419      * corrections.
420      */
421     if (EEL_FULL(fr->eeltype) || EVDW_PME(fr->vdwtype))
422     {
423         real Vlr               = 0, Vcorr = 0;
424         real dvdl_long_range   = 0;
425         int  status            = 0;
426         real Vlr_q             = 0, Vlr_lj = 0, Vcorr_q = 0, Vcorr_lj = 0;
427         real dvdl_long_range_q = 0, dvdl_long_range_lj = 0;
428
429         bSB = (ir->nwall == 2);
430         if (bSB)
431         {
432             copy_mat(box, boxs);
433             svmul(ir->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
434             box_size[ZZ] *= ir->wall_ewald_zfac;
435         }
436
437         if (EEL_PME_EWALD(fr->eeltype) || EVDW_PME(fr->vdwtype))
438         {
439             real dvdl_long_range_correction_q   = 0;
440             real dvdl_long_range_correction_lj  = 0;
441             /* With the Verlet scheme exclusion forces are calculated
442              * in the non-bonded kernel.
443              */
444             /* The TPI molecule does not have exclusions with the rest
445              * of the system and no intra-molecular PME grid
446              * contributions will be calculated in
447              * gmx_pme_calc_energy.
448              */
449             if ((ir->cutoff_scheme == ecutsGROUP && fr->n_tpi == 0) ||
450                 ir->ewald_geometry != eewg3D ||
451                 ir->epsilon_surface != 0)
452             {
453                 int nthreads, t;
454
455                 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
456
457                 if (fr->n_tpi > 0)
458                 {
459                     gmx_fatal(FARGS, "TPI with PME currently only works in a 3D geometry with tin-foil boundary conditions");
460                 }
461
462                 nthreads = gmx_omp_nthreads_get(emntBonded);
463 #pragma omp parallel for num_threads(nthreads) schedule(static)
464                 for (t = 0; t < nthreads; t++)
465                 {
466                     int     s, e, i;
467                     rvec   *fnv;
468                     tensor *vir_q, *vir_lj;
469                     real   *Vcorrt_q, *Vcorrt_lj, *dvdlt_q, *dvdlt_lj;
470                     if (t == 0)
471                     {
472                         fnv       = fr->f_novirsum;
473                         vir_q     = &fr->vir_el_recip;
474                         vir_lj    = &fr->vir_lj_recip;
475                         Vcorrt_q  = &Vcorr_q;
476                         Vcorrt_lj = &Vcorr_lj;
477                         dvdlt_q   = &dvdl_long_range_correction_q;
478                         dvdlt_lj  = &dvdl_long_range_correction_lj;
479                     }
480                     else
481                     {
482                         fnv       = fr->f_t[t].f;
483                         vir_q     = &fr->f_t[t].vir_q;
484                         vir_lj    = &fr->f_t[t].vir_lj;
485                         Vcorrt_q  = &fr->f_t[t].Vcorr_q;
486                         Vcorrt_lj = &fr->f_t[t].Vcorr_lj;
487                         dvdlt_q   = &fr->f_t[t].dvdl[efptCOUL];
488                         dvdlt_lj  = &fr->f_t[t].dvdl[efptVDW];
489                         for (i = 0; i < fr->natoms_force; i++)
490                         {
491                             clear_rvec(fnv[i]);
492                         }
493                         clear_mat(*vir_q);
494                         clear_mat(*vir_lj);
495                     }
496                     *dvdlt_q  = 0;
497                     *dvdlt_lj = 0;
498
499                     ewald_LRcorrection(fr->excl_load[t], fr->excl_load[t+1],
500                                        cr, t, fr,
501                                        md->chargeA,
502                                        md->nChargePerturbed ? md->chargeB : NULL,
503                                        md->sqrt_c6A,
504                                        md->nTypePerturbed ? md->sqrt_c6B : NULL,
505                                        md->sigmaA,
506                                        md->nTypePerturbed ? md->sigmaB : NULL,
507                                        md->sigma3A,
508                                        md->nTypePerturbed ? md->sigma3B : 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_q, *vir_lj,
514                                        Vcorrt_q, Vcorrt_lj,
515                                        lambda[efptCOUL], lambda[efptVDW],
516                                        dvdlt_q, dvdlt_lj);
517                 }
518                 if (nthreads > 1)
519                 {
520                     reduce_thread_forces(fr->natoms_force, fr->f_novirsum,
521                                          fr->vir_el_recip, fr->vir_lj_recip,
522                                          &Vcorr_q, &Vcorr_lj,
523                                          &dvdl_long_range_correction_q,
524                                          &dvdl_long_range_correction_lj,
525                                          nthreads, fr->f_t);
526                 }
527                 wallcycle_sub_stop(wcycle, ewcsEWALD_CORRECTION);
528             }
529
530             if (EEL_PME_EWALD(fr->eeltype) && fr->n_tpi == 0)
531             {
532                 Vcorr_q += ewald_charge_correction(cr, fr, lambda[efptCOUL], box,
533                                                    &dvdl_long_range_correction_q,
534                                                    fr->vir_el_recip);
535             }
536
537             enerd->dvdl_lin[efptCOUL] += dvdl_long_range_correction_q;
538             enerd->dvdl_lin[efptVDW]  += dvdl_long_range_correction_lj;
539
540             if ((EEL_PME(fr->eeltype) || EVDW_PME(fr->vdwtype)) && (cr->duty & DUTY_PME))
541             {
542                 /* Do reciprocal PME for Coulomb and/or LJ. */
543                 assert(fr->n_tpi >= 0);
544                 if (fr->n_tpi == 0 || (flags & GMX_FORCE_STATECHANGED))
545                 {
546                     pme_flags = GMX_PME_SPREAD | GMX_PME_SOLVE;
547                     if (EEL_PME(fr->eeltype))
548                     {
549                         pme_flags     |= GMX_PME_DO_COULOMB;
550                     }
551                     if (EVDW_PME(fr->vdwtype))
552                     {
553                         pme_flags |= GMX_PME_DO_LJ;
554                     }
555                     if (flags & GMX_FORCE_FORCES)
556                     {
557                         pme_flags |= GMX_PME_CALC_F;
558                     }
559                     if (flags & GMX_FORCE_VIRIAL)
560                     {
561                         pme_flags |= GMX_PME_CALC_ENER_VIR;
562                     }
563                     if (fr->n_tpi > 0)
564                     {
565                         /* We don't calculate f, but we do want the potential */
566                         pme_flags |= GMX_PME_CALC_POT;
567                     }
568                     wallcycle_start(wcycle, ewcPMEMESH);
569                     status = gmx_pme_do(fr->pmedata,
570                                         0, md->homenr - fr->n_tpi,
571                                         x, fr->f_novirsum,
572                                         md->chargeA, md->chargeB,
573                                         md->sqrt_c6A, md->sqrt_c6B,
574                                         md->sigmaA, md->sigmaB,
575                                         bSB ? boxs : box, cr,
576                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_x(cr->dd) : 0,
577                                         DOMAINDECOMP(cr) ? dd_pme_maxshift_y(cr->dd) : 0,
578                                         nrnb, wcycle,
579                                         fr->vir_el_recip, fr->ewaldcoeff_q,
580                                         fr->vir_lj_recip, fr->ewaldcoeff_lj,
581                                         &Vlr_q, &Vlr_lj,
582                                         lambda[efptCOUL], lambda[efptVDW],
583                                         &dvdl_long_range_q, &dvdl_long_range_lj, pme_flags);
584                     *cycles_pme = wallcycle_stop(wcycle, ewcPMEMESH);
585                     if (status != 0)
586                     {
587                         gmx_fatal(FARGS, "Error %d in reciprocal PME routine", status);
588                     }
589                     /* We should try to do as little computation after
590                      * this as possible, because parallel PME synchronizes
591                      * the nodes, so we want all load imbalance of the
592                      * rest of the force calculation to be before the PME
593                      * call.  DD load balancing is done on the whole time
594                      * of the force call (without PME).
595                      */
596                 }
597                 if (fr->n_tpi > 0)
598                 {
599                     if (EVDW_PME(ir->vdwtype))
600                     {
601
602                         gmx_fatal(FARGS, "Test particle insertion not implemented with LJ-PME");
603                     }
604                     /* Determine the PME grid energy of the test molecule
605                      * with the PME grid potential of the other charges.
606                      */
607                     gmx_pme_calc_energy(fr->pmedata, fr->n_tpi,
608                                         x + md->homenr - fr->n_tpi,
609                                         md->chargeA + md->homenr - fr->n_tpi,
610                                         &Vlr_q);
611                 }
612             }
613         }
614
615         if (!EEL_PME(fr->eeltype) && EEL_PME_EWALD(fr->eeltype))
616         {
617             Vlr_q = do_ewald(ir, x, fr->f_novirsum,
618                              md->chargeA, md->chargeB,
619                              box_size, cr, md->homenr,
620                              fr->vir_el_recip, fr->ewaldcoeff_q,
621                              lambda[efptCOUL], &dvdl_long_range_q, fr->ewald_table);
622         }
623
624         /* Note that with separate PME nodes we get the real energies later */
625         enerd->dvdl_lin[efptCOUL] += dvdl_long_range_q;
626         enerd->dvdl_lin[efptVDW]  += dvdl_long_range_lj;
627         enerd->term[F_COUL_RECIP]  = Vlr_q + Vcorr_q;
628         enerd->term[F_LJ_RECIP]    = Vlr_lj + Vcorr_lj;
629         if (debug)
630         {
631             fprintf(debug, "Vlr_q = %g, Vcorr_q = %g, Vlr_corr_q = %g\n",
632                     Vlr_q, Vcorr_q, enerd->term[F_COUL_RECIP]);
633             pr_rvecs(debug, 0, "vir_el_recip after corr", fr->vir_el_recip, DIM);
634             pr_rvecs(debug, 0, "fshift after LR Corrections", fr->fshift, SHIFTS);
635             fprintf(debug, "Vlr_lj: %g, Vcorr_lj = %g, Vlr_corr_lj = %g\n",
636                     Vlr_lj, Vcorr_lj, enerd->term[F_LJ_RECIP]);
637             pr_rvecs(debug, 0, "vir_lj_recip after corr", fr->vir_lj_recip, DIM);
638         }
639     }
640     else
641     {
642         /* Is there a reaction-field exclusion correction needed? */
643         if (EEL_RF(fr->eeltype) && eelRF_NEC != fr->eeltype)
644         {
645             /* With the Verlet scheme, exclusion forces are calculated
646              * in the non-bonded kernel.
647              */
648             if (ir->cutoff_scheme != ecutsVERLET)
649             {
650                 real dvdl_rf_excl      = 0;
651                 enerd->term[F_RF_EXCL] =
652                     RF_excl_correction(fr, graph, md, excl, x, f,
653                                        fr->fshift, &pbc, lambda[efptCOUL], &dvdl_rf_excl);
654
655                 enerd->dvdl_lin[efptCOUL] += dvdl_rf_excl;
656             }
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: rank %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 }
692
693 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd)
694 {
695     int i, n2;
696
697     for (i = 0; i < F_NRE; i++)
698     {
699         enerd->term[i]         = 0;
700         enerd->foreign_term[i] = 0;
701     }
702
703
704     for (i = 0; i < efptNR; i++)
705     {
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     enerd->foreign_grpp.nener = n2;
717     for (i = 0; (i < egNR); i++)
718     {
719         snew(enerd->grpp.ener[i], n2);
720         snew(enerd->foreign_grpp.ener[i], n2);
721     }
722
723     if (n_lambda)
724     {
725         enerd->n_lambda = 1 + n_lambda;
726         snew(enerd->enerpart_lambda, enerd->n_lambda);
727     }
728     else
729     {
730         enerd->n_lambda = 0;
731     }
732 }
733
734 void destroy_enerdata(gmx_enerdata_t *enerd)
735 {
736     int i;
737
738     for (i = 0; (i < egNR); i++)
739     {
740         sfree(enerd->grpp.ener[i]);
741     }
742
743     for (i = 0; (i < egNR); i++)
744     {
745         sfree(enerd->foreign_grpp.ener[i]);
746     }
747
748     if (enerd->n_lambda)
749     {
750         sfree(enerd->enerpart_lambda);
751     }
752 }
753
754 static real sum_v(int n, real v[])
755 {
756     real t;
757     int  i;
758
759     t = 0.0;
760     for (i = 0; (i < n); i++)
761     {
762         t = t + v[i];
763     }
764
765     return t;
766 }
767
768 void sum_epot(gmx_grppairener_t *grpp, real *epot)
769 {
770     int i;
771
772     /* Accumulate energies */
773     epot[F_COUL_SR]  = sum_v(grpp->nener, grpp->ener[egCOULSR]);
774     epot[F_LJ]       = sum_v(grpp->nener, grpp->ener[egLJSR]);
775     epot[F_LJ14]     = sum_v(grpp->nener, grpp->ener[egLJ14]);
776     epot[F_COUL14]   = sum_v(grpp->nener, grpp->ener[egCOUL14]);
777     epot[F_COUL_LR]  = sum_v(grpp->nener, grpp->ener[egCOULLR]);
778     epot[F_LJ_LR]    = sum_v(grpp->nener, grpp->ener[egLJLR]);
779     /* We have already added 1-2,1-3, and 1-4 terms to F_GBPOL */
780     epot[F_GBPOL]   += sum_v(grpp->nener, grpp->ener[egGB]);
781
782 /* lattice part of LR doesnt belong to any group
783  * and has been added earlier
784  */
785     epot[F_BHAM]     = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
786     epot[F_BHAM_LR]  = sum_v(grpp->nener, grpp->ener[egBHAMLR]);
787
788     epot[F_EPOT] = 0;
789     for (i = 0; (i < F_EPOT); i++)
790     {
791         if (i != F_DISRESVIOL && i != F_ORIRESDEV)
792         {
793             epot[F_EPOT] += epot[i];
794         }
795     }
796 }
797
798 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals)
799 {
800     int    i, j, index;
801     double dlam;
802
803     enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW];  /* include dispersion correction */
804     enerd->term[F_DVDL]       = 0.0;
805     for (i = 0; i < efptNR; i++)
806     {
807         if (fepvals->separate_dvdl[i])
808         {
809             /* could this be done more readably/compactly? */
810             switch (i)
811             {
812                 case (efptMASS):
813                     index = F_DKDL;
814                     break;
815                 case (efptCOUL):
816                     index = F_DVDL_COUL;
817                     break;
818                 case (efptVDW):
819                     index = F_DVDL_VDW;
820                     break;
821                 case (efptBONDED):
822                     index = F_DVDL_BONDED;
823                     break;
824                 case (efptRESTRAINT):
825                     index = F_DVDL_RESTRAINT;
826                     break;
827                 default:
828                     index = F_DVDL;
829                     break;
830             }
831             enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
832             if (debug)
833             {
834                 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
835                         efpt_names[i], i, enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
836             }
837         }
838         else
839         {
840             enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
841             if (debug)
842             {
843                 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
844                         efpt_names[0], i, enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
845             }
846         }
847     }
848
849     /* Notes on the foreign lambda free energy difference evaluation:
850      * Adding the potential and ekin terms that depend linearly on lambda
851      * as delta lam * dvdl to the energy differences is exact.
852      * For the constraints this is not exact, but we have no other option
853      * without literally changing the lengths and reevaluating the energies at each step.
854      * (try to remedy this post 4.6 - MRS)
855      * For the non-bonded LR term we assume that the soft-core (if present)
856      * no longer affects the energy beyond the short-range cut-off,
857      * which is a very good approximation (except for exotic settings).
858      * (investigate how to overcome this post 4.6 - MRS)
859      */
860     if (fepvals->separate_dvdl[efptBONDED])
861     {
862         enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
863     }
864     else
865     {
866         enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
867     }
868     enerd->term[F_DVDL_CONSTR] = 0;
869
870     for (i = 0; i < fepvals->n_lambda; i++)
871     {
872         /* note we are iterating over fepvals here!
873            For the current lam, dlam = 0 automatically,
874            so we don't need to add anything to the
875            enerd->enerpart_lambda[0] */
876
877         /* we don't need to worry about dvdl_lin contributions to dE at
878            current lambda, because the contributions to the current
879            lambda are automatically zeroed */
880
881         for (j = 0; j < efptNR; j++)
882         {
883             /* Note that this loop is over all dhdl components, not just the separated ones */
884             dlam = (fepvals->all_lambda[j][i]-lambda[j]);
885             enerd->enerpart_lambda[i+1] += dlam*enerd->dvdl_lin[j];
886             if (debug)
887             {
888                 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
889                         fepvals->all_lambda[j][i], efpt_names[j],
890                         (enerd->enerpart_lambda[i+1] - enerd->enerpart_lambda[0]),
891                         dlam, enerd->dvdl_lin[j]);
892             }
893         }
894     }
895 }
896
897
898 void reset_foreign_enerdata(gmx_enerdata_t *enerd)
899 {
900     int  i, j;
901
902     /* First reset all foreign energy components.  Foreign energies always called on
903        neighbor search steps */
904     for (i = 0; (i < egNR); i++)
905     {
906         for (j = 0; (j < enerd->grpp.nener); j++)
907         {
908             enerd->foreign_grpp.ener[i][j] = 0.0;
909         }
910     }
911
912     /* potential energy components */
913     for (i = 0; (i <= F_EPOT); i++)
914     {
915         enerd->foreign_term[i] = 0.0;
916     }
917 }
918
919 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
920                     gmx_enerdata_t *enerd,
921                     gmx_bool bMaster)
922 {
923     gmx_bool bKeepLR;
924     int      i, j;
925
926     /* First reset all energy components, except for the long range terms
927      * on the master at non neighbor search steps, since the long range
928      * terms have already been summed at the last neighbor search step.
929      */
930     bKeepLR = (fr->bTwinRange && !bNS);
931     for (i = 0; (i < egNR); i++)
932     {
933         if (!(bKeepLR && bMaster && (i == egCOULLR || i == egLJLR)))
934         {
935             for (j = 0; (j < enerd->grpp.nener); j++)
936             {
937                 enerd->grpp.ener[i][j] = 0.0;
938             }
939         }
940     }
941     for (i = 0; i < efptNR; i++)
942     {
943         enerd->dvdl_lin[i]    = 0.0;
944         enerd->dvdl_nonlin[i] = 0.0;
945     }
946
947     /* Normal potential energy components */
948     for (i = 0; (i <= F_EPOT); i++)
949     {
950         enerd->term[i] = 0.0;
951     }
952     /* Initialize the dVdlambda term with the long range contribution */
953     /* Initialize the dvdl term with the long range contribution */
954     enerd->term[F_DVDL]            = 0.0;
955     enerd->term[F_DVDL_COUL]       = 0.0;
956     enerd->term[F_DVDL_VDW]        = 0.0;
957     enerd->term[F_DVDL_BONDED]     = 0.0;
958     enerd->term[F_DVDL_RESTRAINT]  = 0.0;
959     enerd->term[F_DKDL]            = 0.0;
960     if (enerd->n_lambda > 0)
961     {
962         for (i = 0; i < enerd->n_lambda; i++)
963         {
964             enerd->enerpart_lambda[i] = 0.0;
965         }
966     }
967     /* reset foreign energy data - separate function since we also call it elsewhere */
968     reset_foreign_enerdata(enerd);
969 }