Replace most strdup with gmx_strdup
[alexxy/gromacs.git] / src / gromacs / mdlib / tpi.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 <math.h>
40 #include <stdlib.h>
41 #include <string.h>
42 #include <time.h>
43
44 #include "network.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "nrnb.h"
47 #include "chargegroup.h"
48 #include "force.h"
49 #include "macros.h"
50 #include "names.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "txtdump.h"
53 #include "typedefs.h"
54 #include "update.h"
55 #include "constr.h"
56 #include "gromacs/math/vec.h"
57 #include "tgroup.h"
58 #include "mdebin.h"
59 #include "vsite.h"
60 #include "force.h"
61 #include "mdrun.h"
62 #include "domdec.h"
63 #include "gromacs/random/random.h"
64 #include "gromacs/math/units.h"
65 #include "gromacs/fileio/xvgr.h"
66 #include "mdatoms.h"
67 #include "ns.h"
68 #include "gromacs/topology/mtop_util.h"
69 #include "pme.h"
70 #include "gromacs/gmxlib/conformation-utilities.h"
71
72 #include "gromacs/legacyheaders/types/commrec.h"
73 #include "gromacs/fileio/confio.h"
74 #include "gromacs/fileio/gmxfio.h"
75 #include "gromacs/fileio/trxio.h"
76 #include "gromacs/timing/wallcycle.h"
77 #include "gromacs/timing/walltime_accounting.h"
78
79 static void global_max(t_commrec *cr, int *n)
80 {
81     int *sum, i;
82
83     snew(sum, cr->nnodes);
84     sum[cr->nodeid] = *n;
85     gmx_sumi(cr->nnodes, sum, cr);
86     for (i = 0; i < cr->nnodes; i++)
87     {
88         *n = max(*n, sum[i]);
89     }
90
91     sfree(sum);
92 }
93
94 static void realloc_bins(double **bin, int *nbin, int nbin_new)
95 {
96     int i;
97
98     if (nbin_new != *nbin)
99     {
100         srenew(*bin, nbin_new);
101         for (i = *nbin; i < nbin_new; i++)
102         {
103             (*bin)[i] = 0;
104         }
105         *nbin = nbin_new;
106     }
107 }
108
109 double do_tpi(FILE *fplog, t_commrec *cr,
110               int nfile, const t_filenm fnm[],
111               const output_env_t oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
112               int gmx_unused nstglobalcomm,
113               gmx_vsite_t gmx_unused *vsite, gmx_constr_t gmx_unused constr,
114               int gmx_unused stepout,
115               t_inputrec *inputrec,
116               gmx_mtop_t *top_global, t_fcdata *fcd,
117               t_state *state,
118               t_mdatoms *mdatoms,
119               t_nrnb *nrnb, gmx_wallcycle_t wcycle,
120               gmx_edsam_t gmx_unused ed,
121               t_forcerec *fr,
122               int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
123               gmx_membed_t gmx_unused membed,
124               real gmx_unused cpt_period, real gmx_unused max_hours,
125               const char gmx_unused *deviceOptions,
126               int gmx_unused imdport,
127               unsigned long gmx_unused Flags,
128               gmx_walltime_accounting_t walltime_accounting)
129 {
130     const char     *TPI = "Test Particle Insertion";
131     gmx_localtop_t *top;
132     gmx_groups_t   *groups;
133     gmx_enerdata_t *enerd;
134     rvec           *f;
135     real            lambda, t, temp, beta, drmax, epot;
136     double          embU, sum_embU, *sum_UgembU, V, V_all, VembU_all;
137     t_trxstatus    *status;
138     t_trxframe      rerun_fr;
139     gmx_bool        bDispCorr, bCharge, bRFExcl, bNotLastFrame, bStateChanged, bNS;
140     tensor          force_vir, shake_vir, vir, pres;
141     int             cg_tp, a_tp0, a_tp1, ngid, gid_tp, nener, e;
142     rvec           *x_mol;
143     rvec            mu_tot, x_init, dx, x_tp;
144     int             nnodes, frame;
145     gmx_int64_t     frame_step_prev, frame_step;
146     gmx_int64_t     nsteps, stepblocksize = 0, step;
147     gmx_int64_t     rnd_count_stride, rnd_count;
148     gmx_int64_t     seed;
149     double          rnd[4];
150     int             i, start, end;
151     FILE           *fp_tpi = NULL;
152     char           *ptr, *dump_pdb, **leg, str[STRLEN], str2[STRLEN];
153     double          dbl, dump_ener;
154     gmx_bool        bCavity;
155     int             nat_cavity  = 0, d;
156     real           *mass_cavity = NULL, mass_tot;
157     int             nbin;
158     double          invbinw, *bin, refvolshift, logV, bUlogV;
159     real            dvdl, prescorr, enercorr, dvdlcorr;
160     gmx_bool        bEnergyOutOfBounds;
161     const char     *tpid_leg[2] = {"direct", "reweighted"};
162
163     /* Since there is no upper limit to the insertion energies,
164      * we need to set an upper limit for the distribution output.
165      */
166     real bU_bin_limit      = 50;
167     real bU_logV_bin_limit = bU_bin_limit + 10;
168
169     nnodes = cr->nnodes;
170
171     top = gmx_mtop_generate_local_top(top_global, inputrec);
172
173     groups = &top_global->groups;
174
175     bCavity = (inputrec->eI == eiTPIC);
176     if (bCavity)
177     {
178         ptr = getenv("GMX_TPIC_MASSES");
179         if (ptr == NULL)
180         {
181             nat_cavity = 1;
182         }
183         else
184         {
185             /* Read (multiple) masses from env var GMX_TPIC_MASSES,
186              * The center of mass of the last atoms is then used for TPIC.
187              */
188             nat_cavity = 0;
189             while (sscanf(ptr, "%lf%n", &dbl, &i) > 0)
190             {
191                 srenew(mass_cavity, nat_cavity+1);
192                 mass_cavity[nat_cavity] = dbl;
193                 fprintf(fplog, "mass[%d] = %f\n",
194                         nat_cavity+1, mass_cavity[nat_cavity]);
195                 nat_cavity++;
196                 ptr += i;
197             }
198             if (nat_cavity == 0)
199             {
200                 gmx_fatal(FARGS, "Found %d masses in GMX_TPIC_MASSES", nat_cavity);
201             }
202         }
203     }
204
205     /*
206        init_em(fplog,TPI,inputrec,&lambda,nrnb,mu_tot,
207        state->box,fr,mdatoms,top,cr,nfile,fnm,NULL,NULL);*/
208     /* We never need full pbc for TPI */
209     fr->ePBC = epbcXYZ;
210     /* Determine the temperature for the Boltzmann weighting */
211     temp = inputrec->opts.ref_t[0];
212     if (fplog)
213     {
214         for (i = 1; (i < inputrec->opts.ngtc); i++)
215         {
216             if (inputrec->opts.ref_t[i] != temp)
217             {
218                 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
219                 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
220             }
221         }
222         fprintf(fplog,
223                 "\n  The temperature for test particle insertion is %.3f K\n\n",
224                 temp);
225     }
226     beta = 1.0/(BOLTZ*temp);
227
228     /* Number of insertions per frame */
229     nsteps = inputrec->nsteps;
230
231     /* Use the same neighborlist with more insertions points
232      * in a sphere of radius drmax around the initial point
233      */
234     /* This should be a proper mdp parameter */
235     drmax = inputrec->rtpi;
236
237     /* An environment variable can be set to dump all configurations
238      * to pdb with an insertion energy <= this value.
239      */
240     dump_pdb  = getenv("GMX_TPI_DUMP");
241     dump_ener = 0;
242     if (dump_pdb)
243     {
244         sscanf(dump_pdb, "%lf", &dump_ener);
245     }
246
247     atoms2md(top_global, inputrec, 0, NULL, top_global->natoms, mdatoms);
248     update_mdatoms(mdatoms, inputrec->fepvals->init_lambda);
249
250     snew(enerd, 1);
251     init_enerdata(groups->grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
252     snew(f, top_global->natoms);
253
254     /* Print to log file  */
255     walltime_accounting_start(walltime_accounting);
256     wallcycle_start(wcycle, ewcRUN);
257     print_start(fplog, cr, walltime_accounting, "Test Particle Insertion");
258
259     /* The last charge group is the group to be inserted */
260     cg_tp = top->cgs.nr - 1;
261     a_tp0 = top->cgs.index[cg_tp];
262     a_tp1 = top->cgs.index[cg_tp+1];
263     if (debug)
264     {
265         fprintf(debug, "TPI cg %d, atoms %d-%d\n", cg_tp, a_tp0, a_tp1);
266     }
267     if (a_tp1 - a_tp0 > 1 &&
268         (inputrec->rlist < inputrec->rcoulomb ||
269          inputrec->rlist < inputrec->rvdw))
270     {
271         gmx_fatal(FARGS, "Can not do TPI for multi-atom molecule with a twin-range cut-off");
272     }
273     snew(x_mol, a_tp1-a_tp0);
274
275     bDispCorr = (inputrec->eDispCorr != edispcNO);
276     bCharge   = FALSE;
277     for (i = a_tp0; i < a_tp1; i++)
278     {
279         /* Copy the coordinates of the molecule to be insterted */
280         copy_rvec(state->x[i], x_mol[i-a_tp0]);
281         /* Check if we need to print electrostatic energies */
282         bCharge |= (mdatoms->chargeA[i] != 0 ||
283                     (mdatoms->chargeB && mdatoms->chargeB[i] != 0));
284     }
285     bRFExcl = (bCharge && EEL_RF(fr->eeltype) && fr->eeltype != eelRF_NEC);
286
287     calc_cgcm(fplog, cg_tp, cg_tp+1, &(top->cgs), state->x, fr->cg_cm);
288     if (bCavity)
289     {
290         if (norm(fr->cg_cm[cg_tp]) > 0.5*inputrec->rlist && fplog)
291         {
292             fprintf(fplog, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
293             fprintf(stderr, "WARNING: Your TPI molecule is not centered at 0,0,0\n");
294         }
295     }
296     else
297     {
298         /* Center the molecule to be inserted at zero */
299         for (i = 0; i < a_tp1-a_tp0; i++)
300         {
301             rvec_dec(x_mol[i], fr->cg_cm[cg_tp]);
302         }
303     }
304
305     if (fplog)
306     {
307         fprintf(fplog, "\nWill insert %d atoms %s partial charges\n",
308                 a_tp1-a_tp0, bCharge ? "with" : "without");
309
310         fprintf(fplog, "\nWill insert %d times in each frame of %s\n",
311                 (int)nsteps, opt2fn("-rerun", nfile, fnm));
312     }
313
314     if (!bCavity)
315     {
316         if (inputrec->nstlist > 1)
317         {
318             if (drmax == 0 && a_tp1-a_tp0 == 1)
319             {
320                 gmx_fatal(FARGS, "Re-using the neighborlist %d times for insertions of a single atom in a sphere of radius %f does not make sense", inputrec->nstlist, drmax);
321             }
322             if (fplog)
323             {
324                 fprintf(fplog, "Will use the same neighborlist for %d insertions in a sphere of radius %f\n", inputrec->nstlist, drmax);
325             }
326         }
327     }
328     else
329     {
330         if (fplog)
331         {
332             fprintf(fplog, "Will insert randomly in a sphere of radius %f around the center of the cavity\n", drmax);
333         }
334     }
335
336     ngid   = groups->grps[egcENER].nr;
337     gid_tp = GET_CGINFO_GID(fr->cginfo[cg_tp]);
338     nener  = 1 + ngid;
339     if (bDispCorr)
340     {
341         nener += 1;
342     }
343     if (bCharge)
344     {
345         nener += ngid;
346         if (bRFExcl)
347         {
348             nener += 1;
349         }
350         if (EEL_FULL(fr->eeltype))
351         {
352             nener += 1;
353         }
354     }
355     snew(sum_UgembU, nener);
356
357     /* Copy the random seed set by the user */
358     seed = inputrec->ld_seed;
359     /* We use the frame step number as one random counter.
360      * The second counter use the insertion (step) count. But we
361      * need multiple random numbers per insertion. This number is
362      * not fixed, since we generate random locations in a sphere
363      * by putting locations in a cube and some of these fail.
364      * A count of 20 is already extremely unlikely, so 10000 is
365      * a safe margin for random numbers per insertion.
366      */
367     rnd_count_stride = 10000;
368
369     if (MASTER(cr))
370     {
371         fp_tpi = xvgropen(opt2fn("-tpi", nfile, fnm),
372                           "TPI energies", "Time (ps)",
373                           "(kJ mol\\S-1\\N) / (nm\\S3\\N)", oenv);
374         xvgr_subtitle(fp_tpi, "f. are averages over one frame", oenv);
375         snew(leg, 4+nener);
376         e = 0;
377         sprintf(str, "-kT log(<Ve\\S-\\betaU\\N>/<V>)");
378         leg[e++] = gmx_strdup(str);
379         sprintf(str, "f. -kT log<e\\S-\\betaU\\N>");
380         leg[e++] = gmx_strdup(str);
381         sprintf(str, "f. <e\\S-\\betaU\\N>");
382         leg[e++] = gmx_strdup(str);
383         sprintf(str, "f. V");
384         leg[e++] = gmx_strdup(str);
385         sprintf(str, "f. <Ue\\S-\\betaU\\N>");
386         leg[e++] = gmx_strdup(str);
387         for (i = 0; i < ngid; i++)
388         {
389             sprintf(str, "f. <U\\sVdW %s\\Ne\\S-\\betaU\\N>",
390                     *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
391             leg[e++] = gmx_strdup(str);
392         }
393         if (bDispCorr)
394         {
395             sprintf(str, "f. <U\\sdisp c\\Ne\\S-\\betaU\\N>");
396             leg[e++] = gmx_strdup(str);
397         }
398         if (bCharge)
399         {
400             for (i = 0; i < ngid; i++)
401             {
402                 sprintf(str, "f. <U\\sCoul %s\\Ne\\S-\\betaU\\N>",
403                         *(groups->grpname[groups->grps[egcENER].nm_ind[i]]));
404                 leg[e++] = gmx_strdup(str);
405             }
406             if (bRFExcl)
407             {
408                 sprintf(str, "f. <U\\sRF excl\\Ne\\S-\\betaU\\N>");
409                 leg[e++] = gmx_strdup(str);
410             }
411             if (EEL_FULL(fr->eeltype))
412             {
413                 sprintf(str, "f. <U\\sCoul recip\\Ne\\S-\\betaU\\N>");
414                 leg[e++] = gmx_strdup(str);
415             }
416         }
417         xvgr_legend(fp_tpi, 4+nener, (const char**)leg, oenv);
418         for (i = 0; i < 4+nener; i++)
419         {
420             sfree(leg[i]);
421         }
422         sfree(leg);
423     }
424     clear_rvec(x_init);
425     V_all     = 0;
426     VembU_all = 0;
427
428     invbinw = 10;
429     nbin    = 10;
430     snew(bin, nbin);
431
432     /* Avoid frame step numbers <= -1 */
433     frame_step_prev = -1;
434
435     bNotLastFrame = read_first_frame(oenv, &status, opt2fn("-rerun", nfile, fnm),
436                                      &rerun_fr, TRX_NEED_X);
437     frame = 0;
438
439     if (rerun_fr.natoms - (bCavity ? nat_cavity : 0) !=
440         mdatoms->nr - (a_tp1 - a_tp0))
441     {
442         gmx_fatal(FARGS, "Number of atoms in trajectory (%d)%s "
443                   "is not equal the number in the run input file (%d) "
444                   "minus the number of atoms to insert (%d)\n",
445                   rerun_fr.natoms, bCavity ? " minus one" : "",
446                   mdatoms->nr, a_tp1-a_tp0);
447     }
448
449     refvolshift = log(det(rerun_fr.box));
450
451     switch (inputrec->eI)
452     {
453         case eiTPI:
454             stepblocksize = inputrec->nstlist;
455             break;
456         case eiTPIC:
457             stepblocksize = 1;
458             break;
459         default:
460             gmx_fatal(FARGS, "Unknown integrator %s", ei_names[inputrec->eI]);
461     }
462
463 #ifdef GMX_SIMD
464     /* Make sure we don't detect SIMD overflow generated before this point */
465     gmx_simd_check_and_reset_overflow();
466 #endif
467
468     while (bNotLastFrame)
469     {
470         frame_step      = rerun_fr.step;
471         if (frame_step <= frame_step_prev)
472         {
473             /* We don't have step number in the trajectory file,
474              * or we have constant or decreasing step numbers.
475              * Ensure we have increasing step numbers, since we use
476              * the step numbers as a counter for random numbers.
477              */
478             frame_step  = frame_step_prev + 1;
479         }
480         frame_step_prev = frame_step;
481
482         lambda = rerun_fr.lambda;
483         t      = rerun_fr.time;
484
485         sum_embU = 0;
486         for (e = 0; e < nener; e++)
487         {
488             sum_UgembU[e] = 0;
489         }
490
491         /* Copy the coordinates from the input trajectory */
492         for (i = 0; i < rerun_fr.natoms; i++)
493         {
494             copy_rvec(rerun_fr.x[i], state->x[i]);
495         }
496         copy_mat(rerun_fr.box, state->box);
497
498         V    = det(state->box);
499         logV = log(V);
500
501         bStateChanged = TRUE;
502         bNS           = TRUE;
503
504         step = cr->nodeid*stepblocksize;
505         while (step < nsteps)
506         {
507             /* Initialize the second counter for random numbers using
508              * the insertion step index. This ensures that we get
509              * the same random numbers independently of how many
510              * MPI ranks we use. Also for the same seed, we get
511              * the same initial random sequence for different nsteps.
512              */
513             rnd_count = step*rnd_count_stride;
514
515             if (!bCavity)
516             {
517                 /* Random insertion in the whole volume */
518                 bNS = (step % inputrec->nstlist == 0);
519                 if (bNS)
520                 {
521                     /* Generate a random position in the box */
522                     gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd);
523                     gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd+2);
524                     for (d = 0; d < DIM; d++)
525                     {
526                         x_init[d] = rnd[d]*state->box[d][d];
527                     }
528                 }
529                 if (inputrec->nstlist == 1)
530                 {
531                     copy_rvec(x_init, x_tp);
532                 }
533                 else
534                 {
535                     /* Generate coordinates within |dx|=drmax of x_init */
536                     do
537                     {
538                         gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd);
539                         gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd+2);
540                         for (d = 0; d < DIM; d++)
541                         {
542                             dx[d] = (2*rnd[d] - 1)*drmax;
543                         }
544                     }
545                     while (norm2(dx) > drmax*drmax);
546                     rvec_add(x_init, dx, x_tp);
547                 }
548             }
549             else
550             {
551                 /* Random insertion around a cavity location
552                  * given by the last coordinate of the trajectory.
553                  */
554                 if (step == 0)
555                 {
556                     if (nat_cavity == 1)
557                     {
558                         /* Copy the location of the cavity */
559                         copy_rvec(rerun_fr.x[rerun_fr.natoms-1], x_init);
560                     }
561                     else
562                     {
563                         /* Determine the center of mass of the last molecule */
564                         clear_rvec(x_init);
565                         mass_tot = 0;
566                         for (i = 0; i < nat_cavity; i++)
567                         {
568                             for (d = 0; d < DIM; d++)
569                             {
570                                 x_init[d] +=
571                                     mass_cavity[i]*rerun_fr.x[rerun_fr.natoms-nat_cavity+i][d];
572                             }
573                             mass_tot += mass_cavity[i];
574                         }
575                         for (d = 0; d < DIM; d++)
576                         {
577                             x_init[d] /= mass_tot;
578                         }
579                     }
580                 }
581                 /* Generate coordinates within |dx|=drmax of x_init */
582                 do
583                 {
584                     gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd);
585                     gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd+2);
586                     for (d = 0; d < DIM; d++)
587                     {
588                         dx[d] = (2*rnd[d] - 1)*drmax;
589                     }
590                 }
591                 while (norm2(dx) > drmax*drmax);
592                 rvec_add(x_init, dx, x_tp);
593             }
594
595             if (a_tp1 - a_tp0 == 1)
596             {
597                 /* Insert a single atom, just copy the insertion location */
598                 copy_rvec(x_tp, state->x[a_tp0]);
599             }
600             else
601             {
602                 /* Copy the coordinates from the top file */
603                 for (i = a_tp0; i < a_tp1; i++)
604                 {
605                     copy_rvec(x_mol[i-a_tp0], state->x[i]);
606                 }
607                 /* Rotate the molecule randomly */
608                 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd);
609                 gmx_rng_cycle_2uniform(frame_step, rnd_count++, seed, RND_SEED_TPI, rnd+2);
610                 rotate_conf(a_tp1-a_tp0, state->x+a_tp0, NULL,
611                             2*M_PI*rnd[0],
612                             2*M_PI*rnd[1],
613                             2*M_PI*rnd[2]);
614                 /* Shift to the insertion location */
615                 for (i = a_tp0; i < a_tp1; i++)
616                 {
617                     rvec_inc(state->x[i], x_tp);
618                 }
619             }
620
621             /* Clear some matrix variables  */
622             clear_mat(force_vir);
623             clear_mat(shake_vir);
624             clear_mat(vir);
625             clear_mat(pres);
626
627             /* Set the charge group center of mass of the test particle */
628             copy_rvec(x_init, fr->cg_cm[top->cgs.nr-1]);
629
630             /* Calc energy (no forces) on new positions.
631              * Since we only need the intermolecular energy
632              * and the RF exclusion terms of the inserted molecule occur
633              * within a single charge group we can pass NULL for the graph.
634              * This also avoids shifts that would move charge groups
635              * out of the box.
636              *
637              * Some checks above ensure than we can not have
638              * twin-range interactions together with nstlist > 1,
639              * therefore we do not need to remember the LR energies.
640              */
641             /* Make do_force do a single node force calculation */
642             cr->nnodes = 1;
643             do_force(fplog, cr, inputrec,
644                      step, nrnb, wcycle, top, &top_global->groups,
645                      state->box, state->x, &state->hist,
646                      f, force_vir, mdatoms, enerd, fcd,
647                      state->lambda,
648                      NULL, fr, NULL, mu_tot, t, NULL, NULL, FALSE,
649                      GMX_FORCE_NONBONDED | GMX_FORCE_ENERGY |
650                      (bNS ? GMX_FORCE_DYNAMICBOX | GMX_FORCE_NS | GMX_FORCE_DO_LR : 0) |
651                      (bStateChanged ? GMX_FORCE_STATECHANGED : 0));
652             cr->nnodes    = nnodes;
653             bStateChanged = FALSE;
654             bNS           = FALSE;
655
656             /* Calculate long range corrections to pressure and energy */
657             calc_dispcorr(fplog, inputrec, fr, step, top_global->natoms, state->box,
658                           lambda, pres, vir, &prescorr, &enercorr, &dvdlcorr);
659             /* figure out how to rearrange the next 4 lines MRS 8/4/2009 */
660             enerd->term[F_DISPCORR]  = enercorr;
661             enerd->term[F_EPOT]     += enercorr;
662             enerd->term[F_PRES]     += prescorr;
663             enerd->term[F_DVDL_VDW] += dvdlcorr;
664
665             epot               = enerd->term[F_EPOT];
666             bEnergyOutOfBounds = FALSE;
667 #ifdef GMX_SIMD_X86_SSE2_OR_HIGHER
668             /* With SSE the energy can overflow, check for this */
669             if (gmx_mm_check_and_reset_overflow())
670             {
671                 if (debug)
672                 {
673                     fprintf(debug, "Found an SSE overflow, assuming the energy is out of bounds\n");
674                 }
675                 bEnergyOutOfBounds = TRUE;
676             }
677 #endif
678             /* If the compiler doesn't optimize this check away
679              * we catch the NAN energies.
680              * The epot>GMX_REAL_MAX check catches inf values,
681              * which should nicely result in embU=0 through the exp below,
682              * but it does not hurt to check anyhow.
683              */
684             /* Non-bonded Interaction usually diverge at r=0.
685              * With tabulated interaction functions the first few entries
686              * should be capped in a consistent fashion between
687              * repulsion, dispersion and Coulomb to avoid accidental
688              * negative values in the total energy.
689              * The table generation code in tables.c does this.
690              * With user tbales the user should take care of this.
691              */
692             if (epot != epot || epot > GMX_REAL_MAX)
693             {
694                 bEnergyOutOfBounds = TRUE;
695             }
696             if (bEnergyOutOfBounds)
697             {
698                 if (debug)
699                 {
700                     fprintf(debug, "\n  time %.3f, step %d: non-finite energy %f, using exp(-bU)=0\n", t, (int)step, epot);
701                 }
702                 embU = 0;
703             }
704             else
705             {
706                 embU      = exp(-beta*epot);
707                 sum_embU += embU;
708                 /* Determine the weighted energy contributions of each energy group */
709                 e                = 0;
710                 sum_UgembU[e++] += epot*embU;
711                 if (fr->bBHAM)
712                 {
713                     for (i = 0; i < ngid; i++)
714                     {
715                         sum_UgembU[e++] +=
716                             (enerd->grpp.ener[egBHAMSR][GID(i, gid_tp, ngid)] +
717                              enerd->grpp.ener[egBHAMLR][GID(i, gid_tp, ngid)])*embU;
718                     }
719                 }
720                 else
721                 {
722                     for (i = 0; i < ngid; i++)
723                     {
724                         sum_UgembU[e++] +=
725                             (enerd->grpp.ener[egLJSR][GID(i, gid_tp, ngid)] +
726                              enerd->grpp.ener[egLJLR][GID(i, gid_tp, ngid)])*embU;
727                     }
728                 }
729                 if (bDispCorr)
730                 {
731                     sum_UgembU[e++] += enerd->term[F_DISPCORR]*embU;
732                 }
733                 if (bCharge)
734                 {
735                     for (i = 0; i < ngid; i++)
736                     {
737                         sum_UgembU[e++] +=
738                             (enerd->grpp.ener[egCOULSR][GID(i, gid_tp, ngid)] +
739                              enerd->grpp.ener[egCOULLR][GID(i, gid_tp, ngid)])*embU;
740                     }
741                     if (bRFExcl)
742                     {
743                         sum_UgembU[e++] += enerd->term[F_RF_EXCL]*embU;
744                     }
745                     if (EEL_FULL(fr->eeltype))
746                     {
747                         sum_UgembU[e++] += enerd->term[F_COUL_RECIP]*embU;
748                     }
749                 }
750             }
751
752             if (embU == 0 || beta*epot > bU_bin_limit)
753             {
754                 bin[0]++;
755             }
756             else
757             {
758                 i = (int)((bU_logV_bin_limit
759                            - (beta*epot - logV + refvolshift))*invbinw
760                           + 0.5);
761                 if (i < 0)
762                 {
763                     i = 0;
764                 }
765                 if (i >= nbin)
766                 {
767                     realloc_bins(&bin, &nbin, i+10);
768                 }
769                 bin[i]++;
770             }
771
772             if (debug)
773             {
774                 fprintf(debug, "TPI %7d %12.5e %12.5f %12.5f %12.5f\n",
775                         (int)step, epot, x_tp[XX], x_tp[YY], x_tp[ZZ]);
776             }
777
778             if (dump_pdb && epot <= dump_ener)
779             {
780                 sprintf(str, "t%g_step%d.pdb", t, (int)step);
781                 sprintf(str2, "t: %f step %d ener: %f", t, (int)step, epot);
782                 write_sto_conf_mtop(str, str2, top_global, state->x, state->v,
783                                     inputrec->ePBC, state->box);
784             }
785
786             step++;
787             if ((step/stepblocksize) % cr->nnodes != cr->nodeid)
788             {
789                 /* Skip all steps assigned to the other MPI ranks */
790                 step += (cr->nnodes - 1)*stepblocksize;
791             }
792         }
793
794         if (PAR(cr))
795         {
796             /* When running in parallel sum the energies over the processes */
797             gmx_sumd(1,    &sum_embU, cr);
798             gmx_sumd(nener, sum_UgembU, cr);
799         }
800
801         frame++;
802         V_all     += V;
803         VembU_all += V*sum_embU/nsteps;
804
805         if (fp_tpi)
806         {
807             if (bVerbose || frame%10 == 0 || frame < 10)
808             {
809                 fprintf(stderr, "mu %10.3e <mu> %10.3e\n",
810                         -log(sum_embU/nsteps)/beta, -log(VembU_all/V_all)/beta);
811             }
812
813             fprintf(fp_tpi, "%10.3f %12.5e %12.5e %12.5e %12.5e",
814                     t,
815                     VembU_all == 0 ? 20/beta : -log(VembU_all/V_all)/beta,
816                     sum_embU == 0  ? 20/beta : -log(sum_embU/nsteps)/beta,
817                     sum_embU/nsteps, V);
818             for (e = 0; e < nener; e++)
819             {
820                 fprintf(fp_tpi, " %12.5e", sum_UgembU[e]/nsteps);
821             }
822             fprintf(fp_tpi, "\n");
823             fflush(fp_tpi);
824         }
825
826         bNotLastFrame = read_next_frame(oenv, status, &rerun_fr);
827     } /* End of the loop  */
828     walltime_accounting_end(walltime_accounting);
829
830     close_trj(status);
831
832     if (fp_tpi != NULL)
833     {
834         gmx_fio_fclose(fp_tpi);
835     }
836
837     if (fplog != NULL)
838     {
839         fprintf(fplog, "\n");
840         fprintf(fplog, "  <V>  = %12.5e nm^3\n", V_all/frame);
841         fprintf(fplog, "  <mu> = %12.5e kJ/mol\n", -log(VembU_all/V_all)/beta);
842     }
843
844     /* Write the Boltzmann factor histogram */
845     if (PAR(cr))
846     {
847         /* When running in parallel sum the bins over the processes */
848         i = nbin;
849         global_max(cr, &i);
850         realloc_bins(&bin, &nbin, i);
851         gmx_sumd(nbin, bin, cr);
852     }
853     if (MASTER(cr))
854     {
855         fp_tpi = xvgropen(opt2fn("-tpid", nfile, fnm),
856                           "TPI energy distribution",
857                           "\\betaU - log(V/<V>)", "count", oenv);
858         sprintf(str, "number \\betaU > %g: %9.3e", bU_bin_limit, bin[0]);
859         xvgr_subtitle(fp_tpi, str, oenv);
860         xvgr_legend(fp_tpi, 2, (const char **)tpid_leg, oenv);
861         for (i = nbin-1; i > 0; i--)
862         {
863             bUlogV = -i/invbinw + bU_logV_bin_limit - refvolshift + log(V_all/frame);
864             fprintf(fp_tpi, "%6.2f %10d %12.5e\n",
865                     bUlogV,
866                     (int)(bin[i]+0.5),
867                     bin[i]*exp(-bUlogV)*V_all/VembU_all);
868         }
869         gmx_fio_fclose(fp_tpi);
870     }
871     sfree(bin);
872
873     sfree(sum_UgembU);
874
875     walltime_accounting_set_nsteps_done(walltime_accounting, frame*inputrec->nsteps);
876
877     return 0;
878 }