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