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