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