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