2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014,2015, 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.
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.
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.
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.
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.
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.
45 #include "thread_mpi/threads.h"
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/domdec/domdec_network.h"
49 #include "gromacs/ewald/pme-load-balancing.h"
50 #include "gromacs/ewald/pme.h"
51 #include "gromacs/fileio/filenm.h"
52 #include "gromacs/fileio/mdoutf.h"
53 #include "gromacs/fileio/trajectory_writing.h"
54 #include "gromacs/fileio/trx.h"
55 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/imd/imd.h"
57 #include "gromacs/legacyheaders/constr.h"
58 #include "gromacs/legacyheaders/ebin.h"
59 #include "gromacs/legacyheaders/force.h"
60 #include "gromacs/legacyheaders/md_logging.h"
61 #include "gromacs/legacyheaders/md_support.h"
62 #include "gromacs/legacyheaders/mdatoms.h"
63 #include "gromacs/legacyheaders/mdebin.h"
64 #include "gromacs/legacyheaders/mdrun.h"
65 #include "gromacs/legacyheaders/network.h"
66 #include "gromacs/legacyheaders/nrnb.h"
67 #include "gromacs/legacyheaders/ns.h"
68 #include "gromacs/legacyheaders/shellfc.h"
69 #include "gromacs/legacyheaders/sighandler.h"
70 #include "gromacs/legacyheaders/sim_util.h"
71 #include "gromacs/legacyheaders/tgroup.h"
72 #include "gromacs/legacyheaders/typedefs.h"
73 #include "gromacs/legacyheaders/update.h"
74 #include "gromacs/legacyheaders/vcm.h"
75 #include "gromacs/legacyheaders/vsite.h"
76 #include "gromacs/legacyheaders/types/commrec.h"
77 #include "gromacs/legacyheaders/types/constr.h"
78 #include "gromacs/legacyheaders/types/enums.h"
79 #include "gromacs/legacyheaders/types/fcdata.h"
80 #include "gromacs/legacyheaders/types/force_flags.h"
81 #include "gromacs/legacyheaders/types/forcerec.h"
82 #include "gromacs/legacyheaders/types/group.h"
83 #include "gromacs/legacyheaders/types/inputrec.h"
84 #include "gromacs/legacyheaders/types/interaction_const.h"
85 #include "gromacs/legacyheaders/types/mdatom.h"
86 #include "gromacs/legacyheaders/types/membedt.h"
87 #include "gromacs/legacyheaders/types/nrnb.h"
88 #include "gromacs/legacyheaders/types/oenv.h"
89 #include "gromacs/legacyheaders/types/shellfc.h"
90 #include "gromacs/legacyheaders/types/state.h"
91 #include "gromacs/listed-forces/manage-threading.h"
92 #include "gromacs/math/utilities.h"
93 #include "gromacs/math/vec.h"
94 #include "gromacs/math/vectypes.h"
95 #include "gromacs/mdlib/compute_io.h"
96 #include "gromacs/mdlib/mdrun_signalling.h"
97 #include "gromacs/mdlib/nb_verlet.h"
98 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
99 #include "gromacs/pbcutil/mshift.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/pulling/pull.h"
102 #include "gromacs/swap/swapcoords.h"
103 #include "gromacs/timing/wallcycle.h"
104 #include "gromacs/timing/walltime_accounting.h"
105 #include "gromacs/topology/atoms.h"
106 #include "gromacs/topology/idef.h"
107 #include "gromacs/topology/mtop_util.h"
108 #include "gromacs/topology/topology.h"
109 #include "gromacs/utility/basedefinitions.h"
110 #include "gromacs/utility/cstringutil.h"
111 #include "gromacs/utility/fatalerror.h"
112 #include "gromacs/utility/real.h"
113 #include "gromacs/utility/smalloc.h"
120 #include "corewrap.h"
123 static void reset_all_counters(FILE *fplog
, t_commrec
*cr
,
125 gmx_int64_t
*step_rel
, t_inputrec
*ir
,
126 gmx_wallcycle_t wcycle
, t_nrnb
*nrnb
,
127 gmx_walltime_accounting_t walltime_accounting
,
128 struct nonbonded_verlet_t
*nbv
)
130 char sbuf
[STEPSTRSIZE
];
132 /* Reset all the counters related to performance over the run */
133 md_print_warn(cr
, fplog
, "step %s: resetting all time and cycle counters\n",
134 gmx_step_str(step
, sbuf
));
138 nbnxn_gpu_reset_timings(nbv
);
141 wallcycle_stop(wcycle
, ewcRUN
);
142 wallcycle_reset_all(wcycle
);
143 if (DOMAINDECOMP(cr
))
145 reset_dd_statistics_counters(cr
->dd
);
148 ir
->init_step
+= *step_rel
;
149 ir
->nsteps
-= *step_rel
;
151 wallcycle_start(wcycle
, ewcRUN
);
152 walltime_accounting_start(walltime_accounting
);
153 print_date_and_time(fplog
, cr
->nodeid
, "Restarted time", gmx_gettime());
156 double do_md(FILE *fplog
, t_commrec
*cr
, int nfile
, const t_filenm fnm
[],
157 const output_env_t oenv
, gmx_bool bVerbose
, gmx_bool bCompact
,
159 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
160 int stepout
, t_inputrec
*ir
,
161 gmx_mtop_t
*top_global
,
163 t_state
*state_global
,
165 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
166 gmx_edsam_t ed
, t_forcerec
*fr
,
167 int repl_ex_nst
, int repl_ex_nex
, int repl_ex_seed
, gmx_membed_t membed
,
168 real cpt_period
, real max_hours
,
171 gmx_walltime_accounting_t walltime_accounting
)
173 gmx_mdoutf_t outf
= NULL
;
174 gmx_int64_t step
, step_rel
;
176 double t
, t0
, lam0
[efptNR
];
177 gmx_bool bGStatEveryStep
, bGStat
, bCalcVir
, bCalcEner
;
178 gmx_bool bNS
, bNStList
, bSimAnn
, bStopCM
, bRerunMD
, bNotLastFrame
= FALSE
,
179 bFirstStep
, bStateFromCP
, bStateFromTPX
, bInitStep
, bLastStep
,
180 bBornRadii
, bStartingFromCpt
;
181 gmx_bool bDoDHDL
= FALSE
, bDoFEP
= FALSE
, bDoExpanded
= FALSE
;
182 gmx_bool do_ene
, do_log
, do_verbose
, bRerunWarnNoV
= TRUE
,
183 bForceUpdate
= FALSE
, bCPT
;
184 gmx_bool bMasterState
;
185 int force_flags
, cglo_flags
;
186 tensor force_vir
, shake_vir
, total_vir
, tmp_vir
, pres
;
193 gmx_repl_ex_t repl_ex
= NULL
;
196 t_mdebin
*mdebin
= NULL
;
197 t_state
*state
= NULL
;
198 rvec
*f_global
= NULL
;
199 gmx_enerdata_t
*enerd
;
201 gmx_global_stat_t gstat
;
202 gmx_update_t upd
= NULL
;
203 t_graph
*graph
= NULL
;
205 gmx_groups_t
*groups
;
206 gmx_ekindata_t
*ekind
;
207 gmx_shellfc_t shellfc
;
208 int count
, nconverged
= 0;
210 gmx_bool bConverged
= TRUE
, bSumEkinhOld
, bDoReplEx
, bExchanged
, bNeedRepartition
;
211 gmx_bool bResetCountersHalfMaxH
= FALSE
;
212 gmx_bool bVV
, bTemp
, bPres
, bTrotter
;
213 gmx_bool bUpdateDoLR
;
222 real saved_conserved_quantity
= 0;
226 char sbuf
[STEPSTRSIZE
], sbuf2
[STEPSTRSIZE
];
227 int handled_stop_condition
= gmx_stop_cond_none
; /* compare to get_stop_condition*/
228 gmx_int64_t multisim_nsteps
= -1; /* number of steps to do before first multisim
229 simulation stops. If equal to zero, don't
230 communicate any more between multisims.*/
231 /* PME load balancing data for GPU kernels */
232 pme_load_balancing_t
*pme_loadbal
= NULL
;
233 gmx_bool bPMETune
= FALSE
;
234 gmx_bool bPMETunePrinting
= FALSE
;
237 gmx_bool bIMDstep
= FALSE
;
240 /* Temporary addition for FAHCORE checkpointing */
244 /* Check for special mdrun options */
245 bRerunMD
= (Flags
& MD_RERUN
);
246 if (Flags
& MD_RESETCOUNTERSHALFWAY
)
250 /* Signal to reset the counters half the simulation steps. */
251 wcycle_set_reset_counters(wcycle
, ir
->nsteps
/2);
253 /* Signal to reset the counters halfway the simulation time. */
254 bResetCountersHalfMaxH
= (max_hours
> 0);
257 /* md-vv uses averaged full step velocities for T-control
258 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
259 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
261 bTrotter
= (bVV
&& (IR_NPT_TROTTER(ir
) || IR_NPH_TROTTER(ir
) || IR_NVT_TROTTER(ir
)));
265 /* Since we don't know if the frames read are related in any way,
266 * rebuild the neighborlist at every step.
269 ir
->nstcalcenergy
= 1;
273 check_ir_old_tpx_versions(cr
, fplog
, ir
, top_global
);
275 nstglobalcomm
= check_nstglobalcomm(fplog
, cr
, nstglobalcomm
, ir
);
276 bGStatEveryStep
= (nstglobalcomm
== 1);
280 ir
->nstxout_compressed
= 0;
282 groups
= &top_global
->groups
;
285 init_md(fplog
, cr
, ir
, oenv
, &t
, &t0
, state_global
->lambda
,
286 &(state_global
->fep_state
), lam0
,
287 nrnb
, top_global
, &upd
,
288 nfile
, fnm
, &outf
, &mdebin
,
289 force_vir
, shake_vir
, mu_tot
, &bSimAnn
, &vcm
, Flags
, wcycle
);
291 clear_mat(total_vir
);
293 /* Energy terms and groups */
295 init_enerdata(top_global
->groups
.grps
[egcENER
].nr
, ir
->fepvals
->n_lambda
,
297 if (DOMAINDECOMP(cr
))
303 snew(f
, top_global
->natoms
);
306 /* Kinetic energy data */
308 init_ekindata(fplog
, top_global
, &(ir
->opts
), ekind
);
309 /* Copy the cos acceleration to the groups struct */
310 ekind
->cosacc
.cos_accel
= ir
->cos_accel
;
312 gstat
= global_stat_init(ir
);
315 /* Check for polarizable models and flexible constraints */
316 shellfc
= init_shell_flexcon(fplog
,
317 top_global
, n_flexible_constraints(constr
),
318 (ir
->bContinuation
||
319 (DOMAINDECOMP(cr
) && !MASTER(cr
))) ?
320 NULL
: state_global
->x
);
321 if (shellfc
&& ir
->nstcalcenergy
!= 1)
323 gmx_fatal(FARGS
, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combinations with shell particles.\nPlease make a new tpr file.", ir
->nstcalcenergy
);
325 if (shellfc
&& DOMAINDECOMP(cr
))
327 gmx_fatal(FARGS
, "Shell particles are not implemented with domain decomposition, use a single rank");
329 if (shellfc
&& ir
->eI
== eiNM
)
331 /* Currently shells don't work with Normal Modes */
332 gmx_fatal(FARGS
, "Normal Mode analysis is not supported with shells.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
335 if (vsite
&& ir
->eI
== eiNM
)
337 /* Currently virtual sites don't work with Normal Modes */
338 gmx_fatal(FARGS
, "Normal Mode analysis is not supported with virtual sites.\nIf you'd like to help with adding support, we have an open discussion at http://redmine.gromacs.org/issues/879\n");
343 tMPI_Thread_mutex_lock(&deform_init_box_mutex
);
344 set_deform_reference_box(upd
,
345 deform_init_init_step_tpx
,
346 deform_init_box_tpx
);
347 tMPI_Thread_mutex_unlock(&deform_init_box_mutex
);
351 double io
= compute_io(ir
, top_global
->natoms
, groups
, mdebin
->ebin
->nener
, 1);
352 if ((io
> 2000) && MASTER(cr
))
355 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
360 if (DOMAINDECOMP(cr
))
362 top
= dd_init_local_top(top_global
);
365 dd_init_local_state(cr
->dd
, state_global
, state
);
367 if (DDMASTER(cr
->dd
) && ir
->nstfout
)
369 snew(f_global
, state_global
->natoms
);
374 top
= gmx_mtop_generate_local_top(top_global
, ir
);
376 forcerec_set_excl_load(fr
, top
);
378 state
= serial_init_local_state(state_global
);
381 atoms2md(top_global
, ir
, 0, NULL
, top_global
->natoms
, mdatoms
);
385 set_vsite_top(vsite
, top
, mdatoms
, cr
);
388 if (ir
->ePBC
!= epbcNONE
&& !fr
->bMolPBC
)
390 graph
= mk_graph(fplog
, &(top
->idef
), 0, top_global
->natoms
, FALSE
, FALSE
);
395 make_local_shells(cr
, mdatoms
, shellfc
);
398 setup_bonded_threading(fr
, &top
->idef
);
401 /* Set up interactive MD (IMD) */
402 init_IMD(ir
, cr
, top_global
, fplog
, ir
->nstcalcenergy
, state_global
->x
,
403 nfile
, fnm
, oenv
, imdport
, Flags
);
405 if (DOMAINDECOMP(cr
))
407 /* Distribute the charge groups over the nodes from the master node */
408 dd_partition_system(fplog
, ir
->init_step
, cr
, TRUE
, 1,
409 state_global
, top_global
, ir
,
410 state
, &f
, mdatoms
, top
, fr
,
411 vsite
, shellfc
, constr
,
416 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
418 if (opt2bSet("-cpi", nfile
, fnm
))
420 bStateFromCP
= gmx_fexist_master(opt2fn_master("-cpi", nfile
, fnm
, cr
), cr
);
424 bStateFromCP
= FALSE
;
429 init_expanded_ensemble(bStateFromCP
, ir
, &state
->dfhist
);
436 /* Update mdebin with energy history if appending to output files */
437 if (Flags
& MD_APPENDFILES
)
439 restore_energyhistory_from_state(mdebin
, &state_global
->enerhist
);
443 /* We might have read an energy history from checkpoint,
444 * free the allocated memory and reset the counts.
446 done_energyhistory(&state_global
->enerhist
);
447 init_energyhistory(&state_global
->enerhist
);
450 /* Set the initial energy history in state by updating once */
451 update_energyhistory(&state_global
->enerhist
, mdebin
);
454 /* Initialize constraints */
455 if (constr
&& !DOMAINDECOMP(cr
))
457 set_constraints(constr
, top
, ir
, mdatoms
, cr
);
460 if (repl_ex_nst
> 0 && MASTER(cr
))
462 repl_ex
= init_replica_exchange(fplog
, cr
->ms
, state_global
, ir
,
463 repl_ex_nst
, repl_ex_nex
, repl_ex_seed
);
466 /* PME tuning is only supported with PME for Coulomb. Is is not supported
467 * with only LJ PME, or for reruns.
469 bPMETune
= ((Flags
& MD_TUNEPME
) && EEL_PME(fr
->eeltype
) && !bRerunMD
&&
470 !(Flags
& MD_REPRODUCIBLE
));
473 pme_loadbal_init(&pme_loadbal
, cr
, fplog
, ir
, state
->box
,
474 fr
->ic
, fr
->pmedata
, use_GPU(fr
->nbv
),
478 if (!ir
->bContinuation
&& !bRerunMD
)
480 if (mdatoms
->cFREEZE
&& (state
->flags
& (1<<estV
)))
482 /* Set the velocities of frozen particles to zero */
483 for (i
= 0; i
< mdatoms
->homenr
; i
++)
485 for (m
= 0; m
< DIM
; m
++)
487 if (ir
->opts
.nFreeze
[mdatoms
->cFREEZE
[i
]][m
])
497 /* Constrain the initial coordinates and velocities */
498 do_constrain_first(fplog
, constr
, ir
, mdatoms
, state
,
503 /* Construct the virtual sites for the initial configuration */
504 construct_vsites(vsite
, state
->x
, ir
->delta_t
, NULL
,
505 top
->idef
.iparams
, top
->idef
.il
,
506 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
512 if (IR_TWINRANGE(*ir
) && repl_ex_nst
% ir
->nstcalclr
!= 0)
514 /* We should exchange at nstcalclr steps to get correct integration */
515 gmx_fatal(FARGS
, "The replica exchange period (%d) is not divisible by nstcalclr (%d)", repl_ex_nst
, ir
->nstcalclr
);
518 if (ir
->efep
!= efepNO
)
520 /* Set free energy calculation frequency as the greatest common
521 * denominator of nstdhdl and repl_ex_nst.
522 * Check for nstcalclr with twin-range, since we need the long-range
523 * contribution to the free-energy at the correct (nstcalclr) steps.
525 nstfep
= ir
->fepvals
->nstdhdl
;
528 if (IR_TWINRANGE(*ir
) &&
529 ir
->expandedvals
->nstexpanded
% ir
->nstcalclr
!= 0)
531 gmx_fatal(FARGS
, "nstexpanded should be divisible by nstcalclr");
533 nstfep
= gmx_greatest_common_divisor(ir
->expandedvals
->nstexpanded
, nstfep
);
537 nstfep
= gmx_greatest_common_divisor(repl_ex_nst
, nstfep
);
539 /* We checked divisibility of repl_ex_nst and nstcalclr above */
540 if (IR_TWINRANGE(*ir
) && nstfep
% ir
->nstcalclr
!= 0)
542 gmx_incons("nstfep not divisible by nstcalclr");
546 /* Be REALLY careful about what flags you set here. You CANNOT assume
547 * this is the first step, since we might be restarting from a checkpoint,
548 * and in that case we should not do any modifications to the state.
550 bStopCM
= (ir
->comm_mode
!= ecmNO
&& !ir
->bContinuation
);
552 cglo_flags
= (CGLO_TEMPERATURE
| CGLO_GSTAT
553 | (bStopCM
? CGLO_STOPCM
: 0)
554 | (bVV
? CGLO_PRESSURE
: 0)
555 | (bVV
? CGLO_CONSTRAINT
: 0)
556 | (bRerunMD
? CGLO_RERUNMD
: 0)
557 | ((Flags
& MD_READ_EKIN
) ? CGLO_READEKIN
: 0));
559 bSumEkinhOld
= FALSE
;
560 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
561 NULL
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
562 constr
, NULL
, FALSE
, state
->box
,
563 top_global
, &bSumEkinhOld
, cglo_flags
);
564 if (ir
->eI
== eiVVAK
)
566 /* a second call to get the half step temperature initialized as well */
567 /* we do the same call as above, but turn the pressure off -- internally to
568 compute_globals, this is recognized as a velocity verlet half-step
569 kinetic energy calculation. This minimized excess variables, but
570 perhaps loses some logic?*/
572 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
573 NULL
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
574 constr
, NULL
, FALSE
, state
->box
,
575 top_global
, &bSumEkinhOld
,
576 cglo_flags
&~(CGLO_STOPCM
| CGLO_PRESSURE
));
579 /* Calculate the initial half step temperature, and save the ekinh_old */
580 if (!(Flags
& MD_STARTFROMCPT
))
582 for (i
= 0; (i
< ir
->opts
.ngtc
); i
++)
584 copy_mat(ekind
->tcstat
[i
].ekinh
, ekind
->tcstat
[i
].ekinh_old
);
589 enerd
->term
[F_TEMP
] *= 2; /* result of averages being done over previous and current step,
590 and there is no previous step */
593 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
594 temperature control */
595 trotter_seq
= init_npt_vars(ir
, state
, &MassQ
, bTrotter
);
599 if (constr
&& !ir
->bContinuation
&& ir
->eConstrAlg
== econtLINCS
)
602 "RMS relative constraint deviation after constraining: %.2e\n",
603 constr_rmsd(constr
, FALSE
));
605 if (EI_STATE_VELOCITY(ir
->eI
))
607 fprintf(fplog
, "Initial temperature: %g K\n", enerd
->term
[F_TEMP
]);
611 fprintf(stderr
, "starting md rerun '%s', reading coordinates from"
612 " input trajectory '%s'\n\n",
613 *(top_global
->name
), opt2fn("-rerun", nfile
, fnm
));
616 fprintf(stderr
, "Calculated time to finish depends on nsteps from "
617 "run input file,\nwhich may not correspond to the time "
618 "needed to process input trajectory.\n\n");
624 fprintf(stderr
, "starting mdrun '%s'\n",
625 *(top_global
->name
));
628 sprintf(tbuf
, "%8.1f", (ir
->init_step
+ir
->nsteps
)*ir
->delta_t
);
632 sprintf(tbuf
, "%s", "infinite");
634 if (ir
->init_step
> 0)
636 fprintf(stderr
, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
637 gmx_step_str(ir
->init_step
+ir
->nsteps
, sbuf
), tbuf
,
638 gmx_step_str(ir
->init_step
, sbuf2
),
639 ir
->init_step
*ir
->delta_t
);
643 fprintf(stderr
, "%s steps, %s ps.\n",
644 gmx_step_str(ir
->nsteps
, sbuf
), tbuf
);
647 fprintf(fplog
, "\n");
650 walltime_accounting_start(walltime_accounting
);
651 wallcycle_start(wcycle
, ewcRUN
);
652 print_start(fplog
, cr
, walltime_accounting
, "mdrun");
654 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
656 chkpt_ret
= fcCheckPointParallel( cr
->nodeid
,
660 gmx_fatal( 3, __FILE__
, __LINE__
, "Checkpoint error on step %d\n", 0 );
665 /***********************************************************
669 ************************************************************/
671 /* if rerunMD then read coordinates and velocities from input trajectory */
674 if (getenv("GMX_FORCE_UPDATE"))
682 bNotLastFrame
= read_first_frame(oenv
, &status
,
683 opt2fn("-rerun", nfile
, fnm
),
684 &rerun_fr
, TRX_NEED_X
| TRX_READ_V
);
685 if (rerun_fr
.natoms
!= top_global
->natoms
)
688 "Number of atoms in trajectory (%d) does not match the "
689 "run input file (%d)\n",
690 rerun_fr
.natoms
, top_global
->natoms
);
692 if (ir
->ePBC
!= epbcNONE
)
696 gmx_fatal(FARGS
, "Rerun trajectory frame step %d time %f does not contain a box, while pbc is used", rerun_fr
.step
, rerun_fr
.time
);
698 if (max_cutoff2(ir
->ePBC
, rerun_fr
.box
) < sqr(fr
->rlistlong
))
700 gmx_fatal(FARGS
, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr
.step
, rerun_fr
.time
);
707 rerun_parallel_comm(cr
, &rerun_fr
, &bNotLastFrame
);
710 if (ir
->ePBC
!= epbcNONE
)
712 /* Set the shift vectors.
713 * Necessary here when have a static box different from the tpr box.
715 calc_shifts(rerun_fr
.box
, fr
->shift_vec
);
719 /* loop over MD steps or if rerunMD to end of input trajectory */
721 /* Skip the first Nose-Hoover integration when we get the state from tpx */
722 bStateFromTPX
= !bStateFromCP
;
723 bInitStep
= bFirstStep
&& (bStateFromTPX
|| bVV
);
724 bStartingFromCpt
= (Flags
& MD_STARTFROMCPT
) && bInitStep
;
725 bSumEkinhOld
= FALSE
;
727 bNeedRepartition
= FALSE
;
729 init_global_signals(&gs
, cr
, ir
, repl_ex_nst
);
731 step
= ir
->init_step
;
734 if (MULTISIM(cr
) && (repl_ex_nst
<= 0 ))
736 /* check how many steps are left in other sims */
737 multisim_nsteps
= get_multisim_nsteps(cr
, ir
->nsteps
);
741 /* and stop now if we should */
742 bLastStep
= (bRerunMD
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
) ||
743 ((multisim_nsteps
>= 0) && (step_rel
>= multisim_nsteps
)));
744 while (!bLastStep
|| (bRerunMD
&& bNotLastFrame
))
747 /* Determine if this is a neighbor search step */
748 bNStList
= (ir
->nstlist
> 0 && step
% ir
->nstlist
== 0);
750 if (bPMETune
&& bNStList
)
752 /* PME grid + cut-off optimization with GPUs or PME nodes */
753 pme_loadbal_do(pme_loadbal
, cr
,
754 (bVerbose
&& MASTER(cr
)) ? stderr
: NULL
,
756 ir
, fr
, state
, wcycle
,
761 wallcycle_start(wcycle
, ewcSTEP
);
767 step
= rerun_fr
.step
;
768 step_rel
= step
- ir
->init_step
;
781 bLastStep
= (step_rel
== ir
->nsteps
);
782 t
= t0
+ step
*ir
->delta_t
;
785 if (ir
->efep
!= efepNO
|| ir
->bSimTemp
)
787 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
788 requiring different logic. */
790 set_current_lambdas(step
, ir
->fepvals
, bRerunMD
, &rerun_fr
, state_global
, state
, lam0
);
791 bDoDHDL
= do_per_step(step
, ir
->fepvals
->nstdhdl
);
792 bDoFEP
= ((ir
->efep
!= efepNO
) && do_per_step(step
, nstfep
));
793 bDoExpanded
= (do_per_step(step
, ir
->expandedvals
->nstexpanded
)
794 && (ir
->bExpanded
) && (step
> 0) && (!bStartingFromCpt
));
797 bDoReplEx
= ((repl_ex_nst
> 0) && (step
> 0) && !bLastStep
&&
798 do_per_step(step
, repl_ex_nst
));
802 update_annealing_target_temp(&(ir
->opts
), t
);
807 if (!DOMAINDECOMP(cr
) || MASTER(cr
))
809 for (i
= 0; i
< state_global
->natoms
; i
++)
811 copy_rvec(rerun_fr
.x
[i
], state_global
->x
[i
]);
815 for (i
= 0; i
< state_global
->natoms
; i
++)
817 copy_rvec(rerun_fr
.v
[i
], state_global
->v
[i
]);
822 for (i
= 0; i
< state_global
->natoms
; i
++)
824 clear_rvec(state_global
->v
[i
]);
828 fprintf(stderr
, "\nWARNING: Some frames do not contain velocities.\n"
829 " Ekin, temperature and pressure are incorrect,\n"
830 " the virial will be incorrect when constraints are present.\n"
832 bRerunWarnNoV
= FALSE
;
836 copy_mat(rerun_fr
.box
, state_global
->box
);
837 copy_mat(state_global
->box
, state
->box
);
839 if (vsite
&& (Flags
& MD_RERUN_VSITE
))
841 if (DOMAINDECOMP(cr
))
843 gmx_fatal(FARGS
, "Vsite recalculation with -rerun is not implemented with domain decomposition, use a single rank");
847 /* Following is necessary because the graph may get out of sync
848 * with the coordinates if we only have every N'th coordinate set
850 mk_mshift(fplog
, graph
, fr
->ePBC
, state
->box
, state
->x
);
851 shift_self(graph
, state
->box
, state
->x
);
853 construct_vsites(vsite
, state
->x
, ir
->delta_t
, state
->v
,
854 top
->idef
.iparams
, top
->idef
.il
,
855 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
858 unshift_self(graph
, state
->box
, state
->x
);
863 /* Stop Center of Mass motion */
864 bStopCM
= (ir
->comm_mode
!= ecmNO
&& do_per_step(step
, ir
->nstcomm
));
868 /* for rerun MD always do Neighbour Searching */
869 bNS
= (bFirstStep
|| ir
->nstlist
!= 0);
874 /* Determine whether or not to do Neighbour Searching and LR */
875 bNS
= (bFirstStep
|| bNStList
|| bExchanged
|| bNeedRepartition
);
878 /* check whether we should stop because another simulation has
882 if ( (multisim_nsteps
>= 0) && (step_rel
>= multisim_nsteps
) &&
883 (multisim_nsteps
!= ir
->nsteps
) )
890 "Stopping simulation %d because another one has finished\n",
894 gs
.sig
[eglsCHKPT
] = 1;
899 /* < 0 means stop at next step, > 0 means stop at next NS step */
900 if ( (gs
.set
[eglsSTOPCOND
] < 0) ||
901 ( (gs
.set
[eglsSTOPCOND
] > 0) && (bNStList
|| ir
->nstlist
== 0) ) )
906 /* Determine whether or not to update the Born radii if doing GB */
907 bBornRadii
= bFirstStep
;
908 if (ir
->implicit_solvent
&& (step
% ir
->nstgbradii
== 0))
913 do_log
= do_per_step(step
, ir
->nstlog
) || bFirstStep
|| bLastStep
;
914 do_verbose
= bVerbose
&&
915 (step
% stepout
== 0 || bFirstStep
|| bLastStep
);
917 if (bNS
&& !(bFirstStep
&& ir
->bContinuation
&& !bRerunMD
))
925 bMasterState
= FALSE
;
926 /* Correct the new box if it is too skewed */
927 if (DYNAMIC_BOX(*ir
))
929 if (correct_box(fplog
, step
, state
->box
, graph
))
934 if (DOMAINDECOMP(cr
) && bMasterState
)
936 dd_collect_state(cr
->dd
, state
, state_global
);
940 if (DOMAINDECOMP(cr
))
942 /* Repartition the domain decomposition */
943 dd_partition_system(fplog
, step
, cr
,
944 bMasterState
, nstglobalcomm
,
945 state_global
, top_global
, ir
,
946 state
, &f
, mdatoms
, top
, fr
,
947 vsite
, shellfc
, constr
,
949 do_verbose
&& !bPMETunePrinting
);
953 if (MASTER(cr
) && do_log
)
955 print_ebin_header(fplog
, step
, t
, state
->lambda
[efptFEP
]); /* can we improve the information printed here? */
958 if (ir
->efep
!= efepNO
)
960 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
963 if ((bRerunMD
&& rerun_fr
.bV
) || bExchanged
)
966 /* We need the kinetic energy at minus the half step for determining
967 * the full step kinetic energy and possibly for T-coupling.*/
968 /* This may not be quite working correctly yet . . . . */
969 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
970 wcycle
, enerd
, NULL
, NULL
, NULL
, NULL
, mu_tot
,
971 constr
, NULL
, FALSE
, state
->box
,
972 top_global
, &bSumEkinhOld
,
973 CGLO_RERUNMD
| CGLO_GSTAT
| CGLO_TEMPERATURE
);
975 clear_mat(force_vir
);
977 /* We write a checkpoint at this MD step when:
978 * either at an NS step when we signalled through gs,
979 * or at the last step (but not when we do not want confout),
980 * but never at the first step or with rerun.
982 bCPT
= (((gs
.set
[eglsCHKPT
] && (bNS
|| ir
->nstlist
== 0)) ||
983 (bLastStep
&& (Flags
& MD_CONFOUT
))) &&
984 step
> ir
->init_step
&& !bRerunMD
);
987 gs
.set
[eglsCHKPT
] = 0;
990 /* Determine the energy and pressure:
991 * at nstcalcenergy steps and at energy output steps (set below).
993 if (EI_VV(ir
->eI
) && (!bInitStep
))
995 /* for vv, the first half of the integration actually corresponds
996 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
997 but the virial needs to be calculated on both the current step and the 'next' step. Future
998 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1000 bCalcEner
= do_per_step(step
-1, ir
->nstcalcenergy
);
1001 bCalcVir
= bCalcEner
||
1002 (ir
->epc
!= epcNO
&& (do_per_step(step
, ir
->nstpcouple
) || do_per_step(step
-1, ir
->nstpcouple
)));
1006 bCalcEner
= do_per_step(step
, ir
->nstcalcenergy
);
1007 bCalcVir
= bCalcEner
||
1008 (ir
->epc
!= epcNO
&& do_per_step(step
, ir
->nstpcouple
));
1011 /* Do we need global communication ? */
1012 bGStat
= (bCalcVir
|| bCalcEner
|| bStopCM
||
1013 do_per_step(step
, nstglobalcomm
) ||
1014 (bVV
&& IR_NVT_TROTTER(ir
) && do_per_step(step
-1, nstglobalcomm
)));
1016 do_ene
= (do_per_step(step
, ir
->nstenergy
) || bLastStep
);
1018 if (do_ene
|| do_log
|| bDoReplEx
)
1025 /* these CGLO_ options remain the same throughout the iteration */
1026 cglo_flags
= ((bRerunMD
? CGLO_RERUNMD
: 0) |
1027 (bGStat
? CGLO_GSTAT
: 0)
1030 force_flags
= (GMX_FORCE_STATECHANGED
|
1031 ((DYNAMIC_BOX(*ir
) || bRerunMD
) ? GMX_FORCE_DYNAMICBOX
: 0) |
1032 GMX_FORCE_ALLFORCES
|
1034 (bCalcVir
? GMX_FORCE_VIRIAL
: 0) |
1035 (bCalcEner
? GMX_FORCE_ENERGY
: 0) |
1036 (bDoFEP
? GMX_FORCE_DHDL
: 0)
1041 if (do_per_step(step
, ir
->nstcalclr
))
1043 force_flags
|= GMX_FORCE_DO_LR
;
1049 /* Now is the time to relax the shells */
1050 count
= relax_shell_flexcon(fplog
, cr
, bVerbose
, step
,
1051 ir
, bNS
, force_flags
,
1054 state
, f
, force_vir
, mdatoms
,
1055 nrnb
, wcycle
, graph
, groups
,
1056 shellfc
, fr
, bBornRadii
, t
, mu_tot
,
1058 mdoutf_get_fp_field(outf
));
1068 /* The coordinates (x) are shifted (to get whole molecules)
1070 * This is parallellized as well, and does communication too.
1071 * Check comments in sim_util.c
1073 do_force(fplog
, cr
, ir
, step
, nrnb
, wcycle
, top
, groups
,
1074 state
->box
, state
->x
, &state
->hist
,
1075 f
, force_vir
, mdatoms
, enerd
, fcd
,
1076 state
->lambda
, graph
,
1077 fr
, vsite
, mu_tot
, t
, mdoutf_get_fp_field(outf
), ed
, bBornRadii
,
1078 (bNS
? GMX_FORCE_NS
: 0) | force_flags
);
1081 if (bVV
&& !bStartingFromCpt
&& !bRerunMD
)
1082 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1086 wallcycle_start(wcycle
, ewcUPDATE
);
1087 if (ir
->eI
== eiVV
&& bInitStep
)
1089 /* if using velocity verlet with full time step Ekin,
1090 * take the first half step only to compute the
1091 * virial for the first step. From there,
1092 * revert back to the initial coordinates
1093 * so that the input is actually the initial step.
1095 snew(vbuf
, state
->natoms
);
1096 copy_rvecn(state
->v
, vbuf
, 0, state
->natoms
); /* should make this better for parallelizing? */
1100 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1101 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ1
);
1104 /* If we are using twin-range interactions where the long-range component
1105 * is only evaluated every nstcalclr>1 steps, we should do a special update
1106 * step to combine the long-range forces on these steps.
1107 * For nstcalclr=1 this is not done, since the forces would have been added
1108 * directly to the short-range forces already.
1110 * TODO Remove various aspects of VV+twin-range in master
1111 * branch, because VV integrators did not ever support
1112 * twin-range multiple time stepping with constraints.
1114 bUpdateDoLR
= (fr
->bTwinRange
&& do_per_step(step
, ir
->nstcalclr
));
1116 update_coords(fplog
, step
, ir
, mdatoms
, state
, fr
->bMolPBC
,
1117 f
, bUpdateDoLR
, fr
->f_twin
, bCalcVir
? &fr
->vir_twin_constr
: NULL
, fcd
,
1118 ekind
, M
, upd
, bInitStep
, etrtVELOCITY1
,
1119 cr
, nrnb
, constr
, &top
->idef
);
1121 if (!bRerunMD
|| rerun_fr
.bV
|| bForceUpdate
) /* Why is rerun_fr.bV here? Unclear. */
1123 wallcycle_stop(wcycle
, ewcUPDATE
);
1124 update_constraints(fplog
, step
, NULL
, ir
, mdatoms
,
1125 state
, fr
->bMolPBC
, graph
, f
,
1126 &top
->idef
, shake_vir
,
1127 cr
, nrnb
, wcycle
, upd
, constr
,
1129 wallcycle_start(wcycle
, ewcUPDATE
);
1130 if (bCalcVir
&& bUpdateDoLR
&& ir
->nstcalclr
> 1)
1132 /* Correct the virial for multiple time stepping */
1133 m_sub(shake_vir
, fr
->vir_twin_constr
, shake_vir
);
1138 /* Need to unshift here if a do_force has been
1139 called in the previous step */
1140 unshift_self(graph
, state
->box
, state
->x
);
1142 /* if VV, compute the pressure and constraints */
1143 /* For VV2, we strictly only need this if using pressure
1144 * control, but we really would like to have accurate pressures
1146 * Think about ways around this in the future?
1147 * For now, keep this choice in comments.
1149 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1150 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1152 bTemp
= ((ir
->eI
== eiVV
&& (!bInitStep
)) || (ir
->eI
== eiVVAK
));
1153 if (bCalcEner
&& ir
->eI
== eiVVAK
)
1155 bSumEkinhOld
= TRUE
;
1157 /* for vv, the first half of the integration actually corresponds to the previous step.
1158 So we need information from the last step in the first half of the integration */
1159 if (bGStat
|| do_per_step(step
-1, nstglobalcomm
))
1161 wallcycle_stop(wcycle
, ewcUPDATE
);
1162 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
1163 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1164 constr
, NULL
, FALSE
, state
->box
,
1165 top_global
, &bSumEkinhOld
,
1168 | (bTemp
? CGLO_TEMPERATURE
: 0)
1169 | (bPres
? CGLO_PRESSURE
: 0)
1170 | (bPres
? CGLO_CONSTRAINT
: 0)
1171 | (bStopCM
? CGLO_STOPCM
: 0)
1174 /* explanation of above:
1175 a) We compute Ekin at the full time step
1176 if 1) we are using the AveVel Ekin, and it's not the
1177 initial step, or 2) if we are using AveEkin, but need the full
1178 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1179 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1180 EkinAveVel because it's needed for the pressure */
1181 wallcycle_start(wcycle
, ewcUPDATE
);
1183 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1188 m_add(force_vir
, shake_vir
, total_vir
); /* we need the un-dispersion corrected total vir here */
1189 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ2
);
1195 wallcycle_stop(wcycle
, ewcUPDATE
);
1196 /* We need the kinetic energy at minus the half step for determining
1197 * the full step kinetic energy and possibly for T-coupling.*/
1198 /* This may not be quite working correctly yet . . . . */
1199 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
1200 wcycle
, enerd
, NULL
, NULL
, NULL
, NULL
, mu_tot
,
1201 constr
, NULL
, FALSE
, state
->box
,
1202 top_global
, &bSumEkinhOld
,
1203 CGLO_RERUNMD
| CGLO_GSTAT
| CGLO_TEMPERATURE
);
1204 wallcycle_start(wcycle
, ewcUPDATE
);
1208 if (bTrotter
&& !bInitStep
)
1210 copy_mat(shake_vir
, state
->svir_prev
);
1211 copy_mat(force_vir
, state
->fvir_prev
);
1212 if (IR_NVT_TROTTER(ir
) && ir
->eI
== eiVV
)
1214 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1215 enerd
->term
[F_TEMP
] = sum_ekin(&(ir
->opts
), ekind
, NULL
, (ir
->eI
== eiVV
), FALSE
);
1216 enerd
->term
[F_EKIN
] = trace(ekind
->ekin
);
1219 /* if it's the initial step, we performed this first step just to get the constraint virial */
1220 if (ir
->eI
== eiVV
&& bInitStep
)
1222 copy_rvecn(vbuf
, state
->v
, 0, state
->natoms
);
1225 wallcycle_stop(wcycle
, ewcUPDATE
);
1228 /* compute the conserved quantity */
1231 saved_conserved_quantity
= compute_conserved_from_auxiliary(ir
, state
, &MassQ
);
1234 last_ekin
= enerd
->term
[F_EKIN
];
1236 if ((ir
->eDispCorr
!= edispcEnerPres
) && (ir
->eDispCorr
!= edispcAllEnerPres
))
1238 saved_conserved_quantity
-= enerd
->term
[F_DISPCORR
];
1240 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1241 if (ir
->efep
!= efepNO
&& !bRerunMD
)
1243 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1247 /* ######## END FIRST UPDATE STEP ############## */
1248 /* ######## If doing VV, we now have v(dt) ###### */
1251 /* perform extended ensemble sampling in lambda - we don't
1252 actually move to the new state before outputting
1253 statistics, but if performing simulated tempering, we
1254 do update the velocities and the tau_t. */
1256 lamnew
= ExpandedEnsembleDynamics(fplog
, ir
, enerd
, state
, &MassQ
, state
->fep_state
, &state
->dfhist
, step
, state
->v
, mdatoms
);
1257 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1258 copy_df_history(&state_global
->dfhist
, &state
->dfhist
);
1261 /* Now we have the energies and forces corresponding to the
1262 * coordinates at time t. We must output all of this before
1265 do_md_trajectory_writing(fplog
, cr
, nfile
, fnm
, step
, step_rel
, t
,
1266 ir
, state
, state_global
, top_global
, fr
,
1267 outf
, mdebin
, ekind
, f
, f_global
,
1269 bCPT
, bRerunMD
, bLastStep
, (Flags
& MD_CONFOUT
),
1271 /* Check if IMD step and do IMD communication, if bIMD is TRUE. */
1272 bIMDstep
= do_IMD(ir
->bIMD
, step
, cr
, bNS
, state
->box
, state
->x
, ir
, t
, wcycle
);
1274 /* kludge -- virial is lost with restart for NPT control. Must restart */
1275 if (bStartingFromCpt
&& bVV
)
1277 copy_mat(state
->svir_prev
, shake_vir
);
1278 copy_mat(state
->fvir_prev
, force_vir
);
1281 elapsed_time
= walltime_accounting_get_current_elapsed_time(walltime_accounting
);
1283 /* Check whether everything is still allright */
1284 if (((int)gmx_get_stop_condition() > handled_stop_condition
)
1285 #ifdef GMX_THREAD_MPI
1290 /* this is just make gs.sig compatible with the hack
1291 of sending signals around by MPI_Reduce with together with
1293 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns
)
1295 gs
.sig
[eglsSTOPCOND
] = 1;
1297 if (gmx_get_stop_condition() == gmx_stop_cond_next
)
1299 gs
.sig
[eglsSTOPCOND
] = -1;
1301 /* < 0 means stop at next step, > 0 means stop at next NS step */
1305 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1306 gmx_get_signal_name(),
1307 gs
.sig
[eglsSTOPCOND
] == 1 ? "NS " : "");
1311 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1312 gmx_get_signal_name(),
1313 gs
.sig
[eglsSTOPCOND
] == 1 ? "NS " : "");
1315 handled_stop_condition
= (int)gmx_get_stop_condition();
1317 else if (MASTER(cr
) && (bNS
|| ir
->nstlist
<= 0) &&
1318 (max_hours
> 0 && elapsed_time
> max_hours
*60.0*60.0*0.99) &&
1319 gs
.sig
[eglsSTOPCOND
] == 0 && gs
.set
[eglsSTOPCOND
] == 0)
1321 /* Signal to terminate the run */
1322 gs
.sig
[eglsSTOPCOND
] = 1;
1325 fprintf(fplog
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step
, sbuf
), max_hours
*0.99);
1327 fprintf(stderr
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step
, sbuf
), max_hours
*0.99);
1330 if (bResetCountersHalfMaxH
&& MASTER(cr
) &&
1331 elapsed_time
> max_hours
*60.0*60.0*0.495)
1333 /* Set flag that will communicate the signal to all ranks in the simulation */
1334 gs
.sig
[eglsRESETCOUNTERS
] = 1;
1337 /* In parallel we only have to check for checkpointing in steps
1338 * where we do global communication,
1339 * otherwise the other nodes don't know.
1341 if (MASTER(cr
) && ((bGStat
|| !PAR(cr
)) &&
1344 elapsed_time
>= nchkpt
*cpt_period
*60.0)) &&
1345 gs
.set
[eglsCHKPT
] == 0)
1347 gs
.sig
[eglsCHKPT
] = 1;
1350 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1355 update_tcouple(step
, ir
, state
, ekind
, &MassQ
, mdatoms
);
1357 if (ETC_ANDERSEN(ir
->etc
)) /* keep this outside of update_tcouple because of the extra info required to pass */
1359 gmx_bool bIfRandomize
;
1360 bIfRandomize
= update_randomize_velocities(ir
, step
, cr
, mdatoms
, state
, upd
, constr
);
1361 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1362 if (constr
&& bIfRandomize
)
1364 update_constraints(fplog
, step
, NULL
, ir
, mdatoms
,
1365 state
, fr
->bMolPBC
, graph
, f
,
1366 &top
->idef
, tmp_vir
,
1367 cr
, nrnb
, wcycle
, upd
, constr
,
1372 /* ######### START SECOND UPDATE STEP ################# */
1373 /* Box is changed in update() when we do pressure coupling,
1374 * but we should still use the old box for energy corrections and when
1375 * writing it to the energy file, so it matches the trajectory files for
1376 * the same timestep above. Make a copy in a separate array.
1378 copy_mat(state
->box
, lastbox
);
1382 if (!bRerunMD
|| rerun_fr
.bV
|| bForceUpdate
)
1384 wallcycle_start(wcycle
, ewcUPDATE
);
1385 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1388 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ3
);
1389 /* We can only do Berendsen coupling after we have summed
1390 * the kinetic energy or virial. Since the happens
1391 * in global_state after update, we should only do it at
1392 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1397 update_tcouple(step
, ir
, state
, ekind
, &MassQ
, mdatoms
);
1398 update_pcouple(fplog
, step
, ir
, state
, pcoupl_mu
, M
, bInitStep
);
1403 bUpdateDoLR
= (fr
->bTwinRange
&& do_per_step(step
, ir
->nstcalclr
));
1405 /* velocity half-step update */
1406 update_coords(fplog
, step
, ir
, mdatoms
, state
, fr
->bMolPBC
, f
,
1407 bUpdateDoLR
, fr
->f_twin
, bCalcVir
? &fr
->vir_twin_constr
: NULL
, fcd
,
1408 ekind
, M
, upd
, FALSE
, etrtVELOCITY2
,
1409 cr
, nrnb
, constr
, &top
->idef
);
1412 /* Above, initialize just copies ekinh into ekin,
1413 * it doesn't copy position (for VV),
1414 * and entire integrator for MD.
1417 if (ir
->eI
== eiVVAK
)
1419 /* We probably only need md->homenr, not state->natoms */
1420 if (state
->natoms
> cbuf_nalloc
)
1422 cbuf_nalloc
= state
->natoms
;
1423 srenew(cbuf
, cbuf_nalloc
);
1425 copy_rvecn(state
->x
, cbuf
, 0, state
->natoms
);
1427 bUpdateDoLR
= (fr
->bTwinRange
&& do_per_step(step
, ir
->nstcalclr
));
1429 update_coords(fplog
, step
, ir
, mdatoms
, state
, fr
->bMolPBC
, f
,
1430 bUpdateDoLR
, fr
->f_twin
, bCalcVir
? &fr
->vir_twin_constr
: NULL
, fcd
,
1431 ekind
, M
, upd
, bInitStep
, etrtPOSITION
, cr
, nrnb
, constr
, &top
->idef
);
1432 wallcycle_stop(wcycle
, ewcUPDATE
);
1434 update_constraints(fplog
, step
, &dvdl_constr
, ir
, mdatoms
, state
,
1435 fr
->bMolPBC
, graph
, f
,
1436 &top
->idef
, shake_vir
,
1437 cr
, nrnb
, wcycle
, upd
, constr
,
1440 if (bCalcVir
&& bUpdateDoLR
&& ir
->nstcalclr
> 1)
1442 /* Correct the virial for multiple time stepping */
1443 m_sub(shake_vir
, fr
->vir_twin_constr
, shake_vir
);
1446 if (ir
->eI
== eiVVAK
)
1448 /* erase F_EKIN and F_TEMP here? */
1449 /* just compute the kinetic energy at the half step to perform a trotter step */
1450 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
1451 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1452 constr
, NULL
, FALSE
, lastbox
,
1453 top_global
, &bSumEkinhOld
,
1454 cglo_flags
| CGLO_TEMPERATURE
1456 wallcycle_start(wcycle
, ewcUPDATE
);
1457 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ4
);
1458 /* now we know the scaling, we can compute the positions again again */
1459 copy_rvecn(cbuf
, state
->x
, 0, state
->natoms
);
1461 bUpdateDoLR
= (fr
->bTwinRange
&& do_per_step(step
, ir
->nstcalclr
));
1463 update_coords(fplog
, step
, ir
, mdatoms
, state
, fr
->bMolPBC
, f
,
1464 bUpdateDoLR
, fr
->f_twin
, bCalcVir
? &fr
->vir_twin_constr
: NULL
, fcd
,
1465 ekind
, M
, upd
, bInitStep
, etrtPOSITION
, cr
, nrnb
, constr
, &top
->idef
);
1466 wallcycle_stop(wcycle
, ewcUPDATE
);
1468 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1469 /* are the small terms in the shake_vir here due
1470 * to numerical errors, or are they important
1471 * physically? I'm thinking they are just errors, but not completely sure.
1472 * For now, will call without actually constraining, constr=NULL*/
1473 update_constraints(fplog
, step
, NULL
, ir
, mdatoms
,
1474 state
, fr
->bMolPBC
, graph
, f
,
1475 &top
->idef
, tmp_vir
,
1476 cr
, nrnb
, wcycle
, upd
, NULL
,
1481 /* this factor or 2 correction is necessary
1482 because half of the constraint force is removed
1483 in the vv step, so we have to double it. See
1484 the Redmine issue #1255. It is not yet clear
1485 if the factor of 2 is exact, or just a very
1486 good approximation, and this will be
1487 investigated. The next step is to see if this
1488 can be done adding a dhdl contribution from the
1489 rattle step, but this is somewhat more
1490 complicated with the current code. Will be
1491 investigated, hopefully for 4.6.3. However,
1492 this current solution is much better than
1493 having it completely wrong.
1495 enerd
->term
[F_DVDL_CONSTR
] += 2*dvdl_constr
;
1499 enerd
->term
[F_DVDL_CONSTR
] += dvdl_constr
;
1504 /* Need to unshift here */
1505 unshift_self(graph
, state
->box
, state
->x
);
1510 wallcycle_start(wcycle
, ewcVSITECONSTR
);
1513 shift_self(graph
, state
->box
, state
->x
);
1515 construct_vsites(vsite
, state
->x
, ir
->delta_t
, state
->v
,
1516 top
->idef
.iparams
, top
->idef
.il
,
1517 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
1521 unshift_self(graph
, state
->box
, state
->x
);
1523 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
1526 /* ############## IF NOT VV, Calculate globals HERE ############ */
1527 /* With Leap-Frog we can skip compute_globals at
1528 * non-communication steps, but we need to calculate
1529 * the kinetic energy one step before communication.
1531 if (bGStat
|| (!EI_VV(ir
->eI
) && do_per_step(step
+1, nstglobalcomm
)))
1533 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
1534 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1536 (step_rel
% gs
.nstms
== 0) &&
1537 (multisim_nsteps
< 0 || (step_rel
< multisim_nsteps
)),
1539 top_global
, &bSumEkinhOld
,
1541 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_ENERGY
: 0)
1542 | (!EI_VV(ir
->eI
) && bStopCM
? CGLO_STOPCM
: 0)
1543 | (!EI_VV(ir
->eI
) ? CGLO_TEMPERATURE
: 0)
1544 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_PRESSURE
: 0)
1549 /* ############# END CALC EKIN AND PRESSURE ################# */
1551 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1552 the virial that should probably be addressed eventually. state->veta has better properies,
1553 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1554 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1556 if (ir
->efep
!= efepNO
&& (!bVV
|| bRerunMD
))
1558 /* Sum up the foreign energy and dhdl terms for md and sd.
1559 Currently done every step so that dhdl is correct in the .edr */
1560 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1562 update_box(fplog
, step
, ir
, mdatoms
, state
, f
,
1563 pcoupl_mu
, nrnb
, upd
);
1565 /* ################# END UPDATE STEP 2 ################# */
1566 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1568 /* The coordinates (x) were unshifted in update */
1571 /* We will not sum ekinh_old,
1572 * so signal that we still have to do it.
1574 bSumEkinhOld
= TRUE
;
1577 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1579 /* use the directly determined last velocity, not actually the averaged half steps */
1580 if (bTrotter
&& ir
->eI
== eiVV
)
1582 enerd
->term
[F_EKIN
] = last_ekin
;
1584 enerd
->term
[F_ETOT
] = enerd
->term
[F_EPOT
] + enerd
->term
[F_EKIN
];
1588 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + saved_conserved_quantity
;
1592 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + compute_conserved_from_auxiliary(ir
, state
, &MassQ
);
1594 /* ######### END PREPARING EDR OUTPUT ########### */
1599 gmx_bool do_dr
, do_or
;
1601 if (fplog
&& do_log
&& bDoExpanded
)
1603 /* only needed if doing expanded ensemble */
1604 PrintFreeEnergyInfoToFile(fplog
, ir
->fepvals
, ir
->expandedvals
, ir
->bSimTemp
? ir
->simtempvals
: NULL
,
1605 &state_global
->dfhist
, state
->fep_state
, ir
->nstlog
, step
);
1607 if (!(bStartingFromCpt
&& (EI_VV(ir
->eI
))))
1611 upd_mdebin(mdebin
, bDoDHDL
, TRUE
,
1612 t
, mdatoms
->tmass
, enerd
, state
,
1613 ir
->fepvals
, ir
->expandedvals
, lastbox
,
1614 shake_vir
, force_vir
, total_vir
, pres
,
1615 ekind
, mu_tot
, constr
);
1619 upd_mdebin_step(mdebin
);
1622 do_dr
= do_per_step(step
, ir
->nstdisreout
);
1623 do_or
= do_per_step(step
, ir
->nstorireout
);
1625 print_ebin(mdoutf_get_fp_ene(outf
), do_ene
, do_dr
, do_or
, do_log
? fplog
: NULL
,
1627 eprNORMAL
, bCompact
, mdebin
, fcd
, groups
, &(ir
->opts
));
1631 pull_print_output(ir
->pull_work
, step
, t
);
1634 if (do_per_step(step
, ir
->nstlog
))
1636 if (fflush(fplog
) != 0)
1638 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
1644 /* Have to do this part _after_ outputting the logfile and the edr file */
1645 /* Gets written into the state at the beginning of next loop*/
1646 state
->fep_state
= lamnew
;
1648 /* Print the remaining wall clock time for the run */
1649 if (MULTIMASTER(cr
) &&
1650 (do_verbose
|| gmx_got_usr_signal()) &&
1655 fprintf(stderr
, "\n");
1657 print_time(stderr
, walltime_accounting
, step
, ir
, cr
);
1660 /* Ion/water position swapping.
1661 * Not done in last step since trajectory writing happens before this call
1662 * in the MD loop and exchanges would be lost anyway. */
1663 bNeedRepartition
= FALSE
;
1664 if ((ir
->eSwapCoords
!= eswapNO
) && (step
> 0) && !bLastStep
&&
1665 do_per_step(step
, ir
->swap
->nstswap
))
1667 bNeedRepartition
= do_swapcoords(cr
, step
, t
, ir
, wcycle
,
1668 bRerunMD
? rerun_fr
.x
: state
->x
,
1669 bRerunMD
? rerun_fr
.box
: state
->box
,
1670 top_global
, MASTER(cr
) && bVerbose
, bRerunMD
);
1672 if (bNeedRepartition
&& DOMAINDECOMP(cr
))
1674 dd_collect_state(cr
->dd
, state
, state_global
);
1678 /* Replica exchange */
1682 bExchanged
= replica_exchange(fplog
, cr
, repl_ex
,
1683 state_global
, enerd
,
1687 if ( (bExchanged
|| bNeedRepartition
) && DOMAINDECOMP(cr
) )
1689 dd_partition_system(fplog
, step
, cr
, TRUE
, 1,
1690 state_global
, top_global
, ir
,
1691 state
, &f
, mdatoms
, top
, fr
,
1692 vsite
, shellfc
, constr
,
1693 nrnb
, wcycle
, FALSE
);
1698 bStartingFromCpt
= FALSE
;
1700 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
1701 /* With all integrators, except VV, we need to retain the pressure
1702 * at the current step for coupling at the next step.
1704 if ((state
->flags
& (1<<estPRES_PREV
)) &&
1706 (ir
->nstpcouple
> 0 && step
% ir
->nstpcouple
== 0)))
1708 /* Store the pressure in t_state for pressure coupling
1709 * at the next MD step.
1711 copy_mat(pres
, state
->pres_prev
);
1714 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
1716 if ( (membed
!= NULL
) && (!bLastStep
) )
1718 rescale_membed(step_rel
, membed
, state_global
->x
);
1725 /* read next frame from input trajectory */
1726 bNotLastFrame
= read_next_frame(oenv
, status
, &rerun_fr
);
1731 rerun_parallel_comm(cr
, &rerun_fr
, &bNotLastFrame
);
1735 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
1736 if (DOMAINDECOMP(cr
) && wcycle
)
1738 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
1741 if (!bRerunMD
|| !rerun_fr
.bStep
)
1743 /* increase the MD step number */
1748 /* TODO make a counter-reset module */
1749 /* If it is time to reset counters, set a flag that remains
1750 true until counters actually get reset */
1751 if (step_rel
== wcycle_get_reset_counters(wcycle
) ||
1752 gs
.set
[eglsRESETCOUNTERS
] != 0)
1754 if (pme_loadbal_is_active(pme_loadbal
))
1756 /* Do not permit counter reset while PME load
1757 * balancing is active. The only purpose for resetting
1758 * counters is to measure reliable performance data,
1759 * and that can't be done before balancing
1762 * TODO consider fixing this by delaying the reset
1763 * until after load balancing completes,
1764 * e.g. https://gerrit.gromacs.org/#/c/4964/2 */
1765 gmx_fatal(FARGS
, "PME tuning was still active when attempting to "
1766 "reset mdrun counters at step %" GMX_PRId64
". Try "
1767 "resetting counters later in the run, e.g. with gmx "
1768 "mdrun -resetstep.", step
);
1770 reset_all_counters(fplog
, cr
, step
, &step_rel
, ir
, wcycle
, nrnb
, walltime_accounting
,
1771 use_GPU(fr
->nbv
) ? fr
->nbv
: NULL
);
1772 wcycle_set_reset_counters(wcycle
, -1);
1773 if (!(cr
->duty
& DUTY_PME
))
1775 /* Tell our PME node to reset its counters */
1776 gmx_pme_send_resetcounters(cr
, step
);
1778 /* Correct max_hours for the elapsed time */
1779 max_hours
-= elapsed_time
/(60.0*60.0);
1780 /* If mdrun -maxh -resethway was active, it can only trigger once */
1781 bResetCountersHalfMaxH
= FALSE
; /* TODO move this to where gs.sig[eglsRESETCOUNTERS] is set */
1782 /* Reset can only happen once, so clear the triggering flag. */
1783 gs
.set
[eglsRESETCOUNTERS
] = 0;
1786 /* If bIMD is TRUE, the master updates the IMD energy record and sends positions to VMD client */
1787 IMD_prep_energies_send_positions(ir
->bIMD
&& MASTER(cr
), bIMDstep
, ir
->imd
, enerd
, step
, bCalcEner
, wcycle
);
1790 /* End of main MD loop */
1793 /* Closing TNG files can include compressing data. Therefore it is good to do that
1794 * before stopping the time measurements. */
1795 mdoutf_tng_close(outf
);
1797 /* Stop measuring walltime */
1798 walltime_accounting_end(walltime_accounting
);
1800 if (bRerunMD
&& MASTER(cr
))
1805 if (!(cr
->duty
& DUTY_PME
))
1807 /* Tell the PME only node to finish */
1808 gmx_pme_send_finish(cr
);
1813 if (ir
->nstcalcenergy
> 0 && !bRerunMD
)
1815 print_ebin(mdoutf_get_fp_ene(outf
), FALSE
, FALSE
, FALSE
, fplog
, step
, t
,
1816 eprAVER
, FALSE
, mdebin
, fcd
, groups
, &(ir
->opts
));
1825 pme_loadbal_done(pme_loadbal
, cr
, fplog
, use_GPU(fr
->nbv
));
1828 if (shellfc
&& fplog
)
1830 fprintf(fplog
, "Fraction of iterations that converged: %.2f %%\n",
1831 (nconverged
*100.0)/step_rel
);
1832 fprintf(fplog
, "Average number of force evaluations per MD step: %.2f\n\n",
1836 if (repl_ex_nst
> 0 && MASTER(cr
))
1838 print_replica_exchange_statistics(fplog
, repl_ex
);
1841 /* IMD cleanup, if bIMD is TRUE. */
1842 IMD_finalize(ir
->bIMD
, ir
->imd
);
1844 walltime_accounting_set_nsteps_done(walltime_accounting
, step_rel
);