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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
58 #include "md_support.h"
59 #include "md_logging.h"
74 #include "mpelogging.h"
76 #include "domdec_network.h"
82 #include "compute_io.h"
84 #include "checkpoint.h"
85 #include "mtop_util.h"
86 #include "sighandler.h"
89 #include "pme_loadbal.h"
92 #include "types/nlistheuristics.h"
93 #include "types/iteratedconstraints.h"
94 #include "nbnxn_cuda_data_mgmt.h"
104 #include "corewrap.h"
107 static void reset_all_counters(FILE *fplog
, t_commrec
*cr
,
108 gmx_large_int_t step
,
109 gmx_large_int_t
*step_rel
, t_inputrec
*ir
,
110 gmx_wallcycle_t wcycle
, t_nrnb
*nrnb
,
111 gmx_runtime_t
*runtime
,
112 nbnxn_cuda_ptr_t cu_nbv
)
114 char sbuf
[STEPSTRSIZE
];
116 /* Reset all the counters related to performance over the run */
117 md_print_warn(cr
, fplog
, "step %s: resetting all time and cycle counters\n",
118 gmx_step_str(step
, sbuf
));
122 nbnxn_cuda_reset_timings(cu_nbv
);
125 wallcycle_stop(wcycle
, ewcRUN
);
126 wallcycle_reset_all(wcycle
);
127 if (DOMAINDECOMP(cr
))
129 reset_dd_statistics_counters(cr
->dd
);
132 ir
->init_step
+= *step_rel
;
133 ir
->nsteps
-= *step_rel
;
135 wallcycle_start(wcycle
, ewcRUN
);
136 runtime_start(runtime
);
137 print_date_and_time(fplog
, cr
->nodeid
, "Restarted time", runtime
);
140 double do_md(FILE *fplog
, t_commrec
*cr
, int nfile
, const t_filenm fnm
[],
141 const output_env_t oenv
, gmx_bool bVerbose
, gmx_bool bCompact
,
143 gmx_vsite_t
*vsite
, gmx_constr_t constr
,
144 int stepout
, t_inputrec
*ir
,
145 gmx_mtop_t
*top_global
,
147 t_state
*state_global
,
149 t_nrnb
*nrnb
, gmx_wallcycle_t wcycle
,
150 gmx_edsam_t ed
, t_forcerec
*fr
,
151 int repl_ex_nst
, int repl_ex_nex
, int repl_ex_seed
, gmx_membed_t membed
,
152 real cpt_period
, real max_hours
,
153 const char *deviceOptions
,
155 gmx_runtime_t
*runtime
)
158 gmx_large_int_t step
, step_rel
;
160 double t
, t0
, lam0
[efptNR
];
161 gmx_bool bGStatEveryStep
, bGStat
, bCalcVir
, bCalcEner
;
162 gmx_bool bNS
, bNStList
, bSimAnn
, bStopCM
, bRerunMD
, bNotLastFrame
= FALSE
,
163 bFirstStep
, bStateFromCP
, bStateFromTPX
, bInitStep
, bLastStep
,
164 bBornRadii
, bStartingFromCpt
;
165 gmx_bool bDoDHDL
= FALSE
, bDoFEP
= FALSE
, bDoExpanded
= FALSE
;
166 gmx_bool do_ene
, do_log
, do_verbose
, bRerunWarnNoV
= TRUE
,
167 bForceUpdate
= FALSE
, bCPT
;
169 gmx_bool bMasterState
;
170 int force_flags
, cglo_flags
;
171 tensor force_vir
, shake_vir
, total_vir
, tmp_vir
, pres
;
176 t_state
*bufstate
= NULL
;
177 matrix
*scale_tot
, pcoupl_mu
, M
, ebox
;
180 gmx_repl_ex_t repl_ex
= NULL
;
183 t_mdebin
*mdebin
= NULL
;
184 t_state
*state
= NULL
;
185 rvec
*f_global
= NULL
;
188 gmx_enerdata_t
*enerd
;
190 gmx_global_stat_t gstat
;
191 gmx_update_t upd
= NULL
;
192 t_graph
*graph
= NULL
;
194 gmx_rng_t mcrng
= NULL
;
196 gmx_groups_t
*groups
;
197 gmx_ekindata_t
*ekind
, *ekind_save
;
198 gmx_shellfc_t shellfc
;
199 int count
, nconverged
= 0;
202 gmx_bool bIonize
= FALSE
;
203 gmx_bool bTCR
= FALSE
, bConverged
= TRUE
, bOK
, bSumEkinhOld
, bExchanged
;
205 gmx_bool bResetCountersHalfMaxH
= FALSE
;
206 gmx_bool bVV
, bIterativeCase
, bFirstIterate
, bTemp
, bPres
, bTrotter
;
207 gmx_bool bUpdateDoLR
;
208 real mu_aver
= 0, dvdl_constr
;
209 int a0
, a1
, gnx
= 0, ii
;
210 atom_id
*grpindex
= NULL
;
212 t_coupl_rec
*tcr
= NULL
;
213 rvec
*xcopy
= NULL
, *vcopy
= NULL
, *cbuf
= NULL
;
214 matrix boxcopy
= {{0}}, lastbox
;
216 real fom
, oldfom
, veta_save
, pcurr
, scalevir
, tracevir
;
223 real saved_conserved_quantity
= 0;
228 char sbuf
[STEPSTRSIZE
], sbuf2
[STEPSTRSIZE
];
229 int handled_stop_condition
= gmx_stop_cond_none
; /* compare to get_stop_condition*/
230 gmx_iterate_t iterate
;
231 gmx_large_int_t multisim_nsteps
= -1; /* number of steps to do before first multisim
232 simulation stops. If equal to zero, don't
233 communicate any more between multisims.*/
234 /* PME load balancing data for GPU kernels */
235 pme_load_balancing_t pme_loadbal
= NULL
;
237 gmx_bool bPMETuneTry
= FALSE
, bPMETuneRunning
= FALSE
;
240 /* Temporary addition for FAHCORE checkpointing */
244 /* Check for special mdrun options */
245 bRerunMD
= (Flags
& MD_RERUN
);
246 bIonize
= (Flags
& MD_IONIZE
);
247 bFFscan
= (Flags
& MD_FFSCAN
);
248 bAppend
= (Flags
& MD_APPENDFILES
);
249 if (Flags
& MD_RESETCOUNTERSHALFWAY
)
253 /* Signal to reset the counters half the simulation steps. */
254 wcycle_set_reset_counters(wcycle
, ir
->nsteps
/2);
256 /* Signal to reset the counters halfway the simulation time. */
257 bResetCountersHalfMaxH
= (max_hours
> 0);
260 /* md-vv uses averaged full step velocities for T-control
261 md-vv-avek uses averaged half step velocities for T-control (but full step ekin for P control)
262 md uses averaged half step kinetic energies to determine temperature unless defined otherwise by GMX_EKIN_AVE_VEL; */
264 if (bVV
) /* to store the initial velocities while computing virial */
266 snew(cbuf
, top_global
->natoms
);
268 /* all the iteratative cases - only if there are constraints */
269 bIterativeCase
= ((IR_NPH_TROTTER(ir
) || IR_NPT_TROTTER(ir
)) && (constr
) && (!bRerunMD
));
270 gmx_iterate_init(&iterate
, FALSE
); /* The default value of iterate->bIterationActive is set to
271 false in this step. The correct value, true or false,
272 is set at each step, as it depends on the frequency of temperature
273 and pressure control.*/
274 bTrotter
= (bVV
&& (IR_NPT_TROTTER(ir
) || IR_NPH_TROTTER(ir
) || IR_NVT_TROTTER(ir
)));
278 /* Since we don't know if the frames read are related in any way,
279 * rebuild the neighborlist at every step.
282 ir
->nstcalcenergy
= 1;
286 check_ir_old_tpx_versions(cr
, fplog
, ir
, top_global
);
288 nstglobalcomm
= check_nstglobalcomm(fplog
, cr
, nstglobalcomm
, ir
);
289 bGStatEveryStep
= (nstglobalcomm
== 1);
291 if (!bGStatEveryStep
&& ir
->nstlist
== -1 && fplog
!= NULL
)
294 "To reduce the energy communication with nstlist = -1\n"
295 "the neighbor list validity should not be checked at every step,\n"
296 "this means that exact integration is not guaranteed.\n"
297 "The neighbor list validity is checked after:\n"
298 " <n.list life time> - 2*std.dev.(n.list life time) steps.\n"
299 "In most cases this will result in exact integration.\n"
300 "This reduces the energy communication by a factor of 2 to 3.\n"
301 "If you want less energy communication, set nstlist > 3.\n\n");
304 if (bRerunMD
|| bFFscan
)
308 groups
= &top_global
->groups
;
311 init_md(fplog
, cr
, ir
, oenv
, &t
, &t0
, state_global
->lambda
,
312 &(state_global
->fep_state
), lam0
,
313 nrnb
, top_global
, &upd
,
314 nfile
, fnm
, &outf
, &mdebin
,
315 force_vir
, shake_vir
, mu_tot
, &bSimAnn
, &vcm
, state_global
, Flags
);
317 clear_mat(total_vir
);
319 /* Energy terms and groups */
321 init_enerdata(top_global
->groups
.grps
[egcENER
].nr
, ir
->fepvals
->n_lambda
,
323 if (DOMAINDECOMP(cr
))
329 snew(f
, top_global
->natoms
);
332 /* Kinetic energy data */
334 init_ekindata(fplog
, top_global
, &(ir
->opts
), ekind
);
335 /* needed for iteration of constraints */
337 init_ekindata(fplog
, top_global
, &(ir
->opts
), ekind_save
);
338 /* Copy the cos acceleration to the groups struct */
339 ekind
->cosacc
.cos_accel
= ir
->cos_accel
;
341 gstat
= global_stat_init(ir
);
344 /* Check for polarizable models and flexible constraints */
345 shellfc
= init_shell_flexcon(fplog
,
346 top_global
, n_flexible_constraints(constr
),
347 (ir
->bContinuation
||
348 (DOMAINDECOMP(cr
) && !MASTER(cr
))) ?
349 NULL
: state_global
->x
);
350 if (shellfc
&& ir
->nstcalcenergy
!= 1)
352 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
);
354 if (shellfc
&& DOMAINDECOMP(cr
))
356 gmx_fatal(FARGS
, "In order to run parallel simulations with shells you need to use the -pd flag to mdrun.");
358 if (shellfc
&& ir
->eI
== eiNM
)
360 /* Currently shells don't work with Normal Modes */
361 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");
364 if (vsite
&& ir
->eI
== eiNM
)
366 /* Currently virtual sites don't work with Normal Modes */
367 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");
372 #ifdef GMX_THREAD_MPI
373 tMPI_Thread_mutex_lock(&deform_init_box_mutex
);
375 set_deform_reference_box(upd
,
376 deform_init_init_step_tpx
,
377 deform_init_box_tpx
);
378 #ifdef GMX_THREAD_MPI
379 tMPI_Thread_mutex_unlock(&deform_init_box_mutex
);
384 double io
= compute_io(ir
, top_global
->natoms
, groups
, mdebin
->ebin
->nener
, 1);
385 if ((io
> 2000) && MASTER(cr
))
388 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
393 if (DOMAINDECOMP(cr
))
395 top
= dd_init_local_top(top_global
);
398 dd_init_local_state(cr
->dd
, state_global
, state
);
400 if (DDMASTER(cr
->dd
) && ir
->nstfout
)
402 snew(f_global
, state_global
->natoms
);
409 /* Initialize the particle decomposition and split the topology */
410 top
= split_system(fplog
, top_global
, ir
, cr
);
412 pd_cg_range(cr
, &fr
->cg0
, &fr
->hcg
);
413 pd_at_range(cr
, &a0
, &a1
);
417 top
= gmx_mtop_generate_local_top(top_global
, ir
);
420 a1
= top_global
->natoms
;
423 forcerec_set_excl_load(fr
, top
, cr
);
425 state
= partdec_init_local_state(cr
, state_global
);
428 atoms2md(top_global
, ir
, 0, NULL
, a0
, a1
-a0
, mdatoms
);
432 set_vsite_top(vsite
, top
, mdatoms
, cr
);
435 if (ir
->ePBC
!= epbcNONE
&& !fr
->bMolPBC
)
437 graph
= mk_graph(fplog
, &(top
->idef
), 0, top_global
->natoms
, FALSE
, FALSE
);
442 make_local_shells(cr
, mdatoms
, shellfc
);
445 setup_bonded_threading(fr
, &top
->idef
);
447 if (ir
->pull
&& PAR(cr
))
449 dd_make_local_pull_groups(NULL
, ir
->pull
, mdatoms
);
453 if (DOMAINDECOMP(cr
))
455 /* Distribute the charge groups over the nodes from the master node */
456 dd_partition_system(fplog
, ir
->init_step
, cr
, TRUE
, 1,
457 state_global
, top_global
, ir
,
458 state
, &f
, mdatoms
, top
, fr
,
459 vsite
, shellfc
, constr
,
460 nrnb
, wcycle
, FALSE
);
464 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
466 if (opt2bSet("-cpi", nfile
, fnm
))
468 bStateFromCP
= gmx_fexist_master(opt2fn_master("-cpi", nfile
, fnm
, cr
), cr
);
472 bStateFromCP
= FALSE
;
477 init_expanded_ensemble(bStateFromCP
,ir
,&mcrng
,&state
->dfhist
);
484 /* Update mdebin with energy history if appending to output files */
485 if (Flags
& MD_APPENDFILES
)
487 restore_energyhistory_from_state(mdebin
, &state_global
->enerhist
);
491 /* We might have read an energy history from checkpoint,
492 * free the allocated memory and reset the counts.
494 done_energyhistory(&state_global
->enerhist
);
495 init_energyhistory(&state_global
->enerhist
);
498 /* Set the initial energy history in state by updating once */
499 update_energyhistory(&state_global
->enerhist
, mdebin
);
502 if ((state
->flags
& (1<<estLD_RNG
)) && (Flags
& MD_READ_RNG
))
504 /* Set the random state if we read a checkpoint file */
505 set_stochd_state(upd
, state
);
508 if (state
->flags
& (1<<estMC_RNG
))
510 set_mc_state(mcrng
, state
);
513 /* Initialize constraints */
516 if (!DOMAINDECOMP(cr
))
518 set_constraints(constr
, top
, ir
, mdatoms
, cr
);
522 /* Check whether we have to GCT stuff */
523 bTCR
= ftp2bSet(efGCT
, nfile
, fnm
);
528 fprintf(stderr
, "Will do General Coupling Theory!\n");
530 gnx
= top_global
->mols
.nr
;
532 for (i
= 0; (i
< gnx
); i
++)
540 /* We need to be sure replica exchange can only occur
541 * when the energies are current */
542 check_nst_param(fplog
, cr
, "nstcalcenergy", ir
->nstcalcenergy
,
543 "repl_ex_nst", &repl_ex_nst
);
544 /* This check needs to happen before inter-simulation
545 * signals are initialized, too */
547 if (repl_ex_nst
> 0 && MASTER(cr
))
549 repl_ex
= init_replica_exchange(fplog
, cr
->ms
, state_global
, ir
,
550 repl_ex_nst
, repl_ex_nex
, repl_ex_seed
);
553 /* PME tuning is only supported with GPUs or PME nodes and not with rerun.
554 * With perturbed charges with soft-core we should not change the cut-off.
556 if ((Flags
& MD_TUNEPME
) &&
557 EEL_PME(fr
->eeltype
) &&
558 ( (fr
->cutoff_scheme
== ecutsVERLET
&& fr
->nbv
->bUseGPU
) || !(cr
->duty
& DUTY_PME
)) &&
559 !(ir
->efep
!= efepNO
&& mdatoms
->nChargePerturbed
> 0 && ir
->fepvals
->bScCoul
) &&
562 pme_loadbal_init(&pme_loadbal
, ir
, state
->box
, fr
->ic
, fr
->pmedata
);
564 if (cr
->duty
& DUTY_PME
)
566 /* Start tuning right away, as we can't measure the load */
567 bPMETuneRunning
= TRUE
;
571 /* Separate PME nodes, we can measure the PP/PME load balance */
576 if (!ir
->bContinuation
&& !bRerunMD
)
578 if (mdatoms
->cFREEZE
&& (state
->flags
& (1<<estV
)))
580 /* Set the velocities of frozen particles to zero */
581 for (i
= mdatoms
->start
; i
< mdatoms
->start
+mdatoms
->homenr
; i
++)
583 for (m
= 0; m
< DIM
; m
++)
585 if (ir
->opts
.nFreeze
[mdatoms
->cFREEZE
[i
]][m
])
595 /* Constrain the initial coordinates and velocities */
596 do_constrain_first(fplog
, constr
, ir
, mdatoms
, state
, f
,
597 graph
, cr
, nrnb
, fr
, top
, shake_vir
);
601 /* Construct the virtual sites for the initial configuration */
602 construct_vsites(fplog
, vsite
, state
->x
, nrnb
, ir
->delta_t
, NULL
,
603 top
->idef
.iparams
, top
->idef
.il
,
604 fr
->ePBC
, fr
->bMolPBC
, graph
, cr
, state
->box
);
610 /* set free energy calculation frequency as the minimum
611 greatest common denominator of nstdhdl, nstexpanded, and repl_ex_nst*/
612 nstfep
= ir
->fepvals
->nstdhdl
;
615 nstfep
= gmx_greatest_common_divisor(ir
->fepvals
->nstdhdl
,nstfep
);
619 nstfep
= gmx_greatest_common_divisor(repl_ex_nst
,nstfep
);
622 /* I'm assuming we need global communication the first time! MRS */
623 cglo_flags
= (CGLO_TEMPERATURE
| CGLO_GSTAT
624 | ((ir
->comm_mode
!= ecmNO
) ? CGLO_STOPCM
: 0)
625 | (bVV
? CGLO_PRESSURE
: 0)
626 | (bVV
? CGLO_CONSTRAINT
: 0)
627 | (bRerunMD
? CGLO_RERUNMD
: 0)
628 | ((Flags
& MD_READ_EKIN
) ? CGLO_READEKIN
: 0));
630 bSumEkinhOld
= FALSE
;
631 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
632 NULL
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
633 constr
, NULL
, FALSE
, state
->box
,
634 top_global
, &pcurr
, top_global
->natoms
, &bSumEkinhOld
, cglo_flags
);
635 if (ir
->eI
== eiVVAK
)
637 /* a second call to get the half step temperature initialized as well */
638 /* we do the same call as above, but turn the pressure off -- internally to
639 compute_globals, this is recognized as a velocity verlet half-step
640 kinetic energy calculation. This minimized excess variables, but
641 perhaps loses some logic?*/
643 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
644 NULL
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
645 constr
, NULL
, FALSE
, state
->box
,
646 top_global
, &pcurr
, top_global
->natoms
, &bSumEkinhOld
,
647 cglo_flags
&~(CGLO_STOPCM
| CGLO_PRESSURE
));
650 /* Calculate the initial half step temperature, and save the ekinh_old */
651 if (!(Flags
& MD_STARTFROMCPT
))
653 for (i
= 0; (i
< ir
->opts
.ngtc
); i
++)
655 copy_mat(ekind
->tcstat
[i
].ekinh
, ekind
->tcstat
[i
].ekinh_old
);
660 enerd
->term
[F_TEMP
] *= 2; /* result of averages being done over previous and current step,
661 and there is no previous step */
664 /* if using an iterative algorithm, we need to create a working directory for the state. */
667 bufstate
= init_bufstate(state
);
671 snew(xcopy
, state
->natoms
);
672 snew(vcopy
, state
->natoms
);
673 copy_rvecn(state
->x
, xcopy
, 0, state
->natoms
);
674 copy_rvecn(state
->v
, vcopy
, 0, state
->natoms
);
675 copy_mat(state
->box
, boxcopy
);
678 /* need to make an initiation call to get the Trotter variables set, as well as other constants for non-trotter
679 temperature control */
680 trotter_seq
= init_npt_vars(ir
, state
, &MassQ
, bTrotter
);
684 if (constr
&& !ir
->bContinuation
&& ir
->eConstrAlg
== econtLINCS
)
687 "RMS relative constraint deviation after constraining: %.2e\n",
688 constr_rmsd(constr
, FALSE
));
690 if (EI_STATE_VELOCITY(ir
->eI
))
692 fprintf(fplog
, "Initial temperature: %g K\n", enerd
->term
[F_TEMP
]);
696 fprintf(stderr
, "starting md rerun '%s', reading coordinates from"
697 " input trajectory '%s'\n\n",
698 *(top_global
->name
), opt2fn("-rerun", nfile
, fnm
));
701 fprintf(stderr
, "Calculated time to finish depends on nsteps from "
702 "run input file,\nwhich may not correspond to the time "
703 "needed to process input trajectory.\n\n");
709 fprintf(stderr
, "starting mdrun '%s'\n",
710 *(top_global
->name
));
713 sprintf(tbuf
, "%8.1f", (ir
->init_step
+ir
->nsteps
)*ir
->delta_t
);
717 sprintf(tbuf
, "%s", "infinite");
719 if (ir
->init_step
> 0)
721 fprintf(stderr
, "%s steps, %s ps (continuing from step %s, %8.1f ps).\n",
722 gmx_step_str(ir
->init_step
+ir
->nsteps
, sbuf
), tbuf
,
723 gmx_step_str(ir
->init_step
, sbuf2
),
724 ir
->init_step
*ir
->delta_t
);
728 fprintf(stderr
, "%s steps, %s ps.\n",
729 gmx_step_str(ir
->nsteps
, sbuf
), tbuf
);
732 fprintf(fplog
, "\n");
735 print_start(fplog
, cr
, runtime
, "mdrun");
736 runtime_start(runtime
);
737 wallcycle_start(wcycle
, ewcRUN
);
739 /* safest point to do file checkpointing is here. More general point would be immediately before integrator call */
741 chkpt_ret
= fcCheckPointParallel( cr
->nodeid
,
745 gmx_fatal( 3, __FILE__
, __LINE__
, "Checkpoint error on step %d\n", 0 );
750 /***********************************************************
754 ************************************************************/
756 /* if rerunMD then read coordinates and velocities from input trajectory */
759 if (getenv("GMX_FORCE_UPDATE"))
767 bNotLastFrame
= read_first_frame(oenv
, &status
,
768 opt2fn("-rerun", nfile
, fnm
),
769 &rerun_fr
, TRX_NEED_X
| TRX_READ_V
);
770 if (rerun_fr
.natoms
!= top_global
->natoms
)
773 "Number of atoms in trajectory (%d) does not match the "
774 "run input file (%d)\n",
775 rerun_fr
.natoms
, top_global
->natoms
);
777 if (ir
->ePBC
!= epbcNONE
)
781 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
);
783 if (max_cutoff2(ir
->ePBC
, rerun_fr
.box
) < sqr(fr
->rlistlong
))
785 gmx_fatal(FARGS
, "Rerun trajectory frame step %d time %f has too small box dimensions", rerun_fr
.step
, rerun_fr
.time
);
792 rerun_parallel_comm(cr
, &rerun_fr
, &bNotLastFrame
);
795 if (ir
->ePBC
!= epbcNONE
)
797 /* Set the shift vectors.
798 * Necessary here when have a static box different from the tpr box.
800 calc_shifts(rerun_fr
.box
, fr
->shift_vec
);
804 /* loop over MD steps or if rerunMD to end of input trajectory */
806 /* Skip the first Nose-Hoover integration when we get the state from tpx */
807 bStateFromTPX
= !bStateFromCP
;
808 bInitStep
= bFirstStep
&& (bStateFromTPX
|| bVV
);
809 bStartingFromCpt
= (Flags
& MD_STARTFROMCPT
) && bInitStep
;
811 bSumEkinhOld
= FALSE
;
814 init_global_signals(&gs
, cr
, ir
, repl_ex_nst
);
816 step
= ir
->init_step
;
819 if (ir
->nstlist
== -1)
821 init_nlistheuristics(&nlh
, bGStatEveryStep
, step
);
824 if (MULTISIM(cr
) && (repl_ex_nst
<= 0 ))
826 /* check how many steps are left in other sims */
827 multisim_nsteps
= get_multisim_nsteps(cr
, ir
->nsteps
);
831 /* and stop now if we should */
832 bLastStep
= (bRerunMD
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
) ||
833 ((multisim_nsteps
>= 0) && (step_rel
>= multisim_nsteps
)));
834 while (!bLastStep
|| (bRerunMD
&& bNotLastFrame
))
837 wallcycle_start(wcycle
, ewcSTEP
);
839 GMX_MPE_LOG(ev_timestep1
);
845 step
= rerun_fr
.step
;
846 step_rel
= step
- ir
->init_step
;
859 bLastStep
= (step_rel
== ir
->nsteps
);
860 t
= t0
+ step
*ir
->delta_t
;
863 if (ir
->efep
!= efepNO
|| ir
->bSimTemp
)
865 /* find and set the current lambdas. If rerunning, we either read in a state, or a lambda value,
866 requiring different logic. */
868 set_current_lambdas(step
, ir
->fepvals
, bRerunMD
, &rerun_fr
, state_global
, state
, lam0
);
869 bDoDHDL
= do_per_step(step
, ir
->fepvals
->nstdhdl
);
870 bDoFEP
= (do_per_step(step
, nstfep
) && (ir
->efep
!= efepNO
));
871 bDoExpanded
= (do_per_step(step
, ir
->expandedvals
->nstexpanded
)
872 && (ir
->bExpanded
) && (step
> 0) && (!bStartingFromCpt
));
877 update_annealing_target_temp(&(ir
->opts
), t
);
882 if (!(DOMAINDECOMP(cr
) && !MASTER(cr
)))
884 for (i
= 0; i
< state_global
->natoms
; i
++)
886 copy_rvec(rerun_fr
.x
[i
], state_global
->x
[i
]);
890 for (i
= 0; i
< state_global
->natoms
; i
++)
892 copy_rvec(rerun_fr
.v
[i
], state_global
->v
[i
]);
897 for (i
= 0; i
< state_global
->natoms
; i
++)
899 clear_rvec(state_global
->v
[i
]);
903 fprintf(stderr
, "\nWARNING: Some frames do not contain velocities.\n"
904 " Ekin, temperature and pressure are incorrect,\n"
905 " the virial will be incorrect when constraints are present.\n"
907 bRerunWarnNoV
= FALSE
;
911 copy_mat(rerun_fr
.box
, state_global
->box
);
912 copy_mat(state_global
->box
, state
->box
);
914 if (vsite
&& (Flags
& MD_RERUN_VSITE
))
916 if (DOMAINDECOMP(cr
))
918 gmx_fatal(FARGS
, "Vsite recalculation with -rerun is not implemented for domain decomposition, use particle decomposition");
922 /* Following is necessary because the graph may get out of sync
923 * with the coordinates if we only have every N'th coordinate set
925 mk_mshift(fplog
, graph
, fr
->ePBC
, state
->box
, state
->x
);
926 shift_self(graph
, state
->box
, state
->x
);
928 construct_vsites(fplog
, vsite
, state
->x
, nrnb
, ir
->delta_t
, state
->v
,
929 top
->idef
.iparams
, top
->idef
.il
,
930 fr
->ePBC
, fr
->bMolPBC
, graph
, cr
, state
->box
);
933 unshift_self(graph
, state
->box
, state
->x
);
938 /* Stop Center of Mass motion */
939 bStopCM
= (ir
->comm_mode
!= ecmNO
&& do_per_step(step
, ir
->nstcomm
));
941 /* Copy back starting coordinates in case we're doing a forcefield scan */
944 for (ii
= 0; (ii
< state
->natoms
); ii
++)
946 copy_rvec(xcopy
[ii
], state
->x
[ii
]);
947 copy_rvec(vcopy
[ii
], state
->v
[ii
]);
949 copy_mat(boxcopy
, state
->box
);
954 /* for rerun MD always do Neighbour Searching */
955 bNS
= (bFirstStep
|| ir
->nstlist
!= 0);
960 /* Determine whether or not to do Neighbour Searching and LR */
961 bNStList
= (ir
->nstlist
> 0 && step
% ir
->nstlist
== 0);
963 bNS
= (bFirstStep
|| bExchanged
|| bNStList
|| bDoFEP
||
964 (ir
->nstlist
== -1 && nlh
.nabnsb
> 0));
966 if (bNS
&& ir
->nstlist
== -1)
968 set_nlistheuristics(&nlh
, bFirstStep
|| bExchanged
|| bDoFEP
, step
);
972 /* check whether we should stop because another simulation has
976 if ( (multisim_nsteps
>= 0) && (step_rel
>= multisim_nsteps
) &&
977 (multisim_nsteps
!= ir
->nsteps
) )
984 "Stopping simulation %d because another one has finished\n",
988 gs
.sig
[eglsCHKPT
] = 1;
993 /* < 0 means stop at next step, > 0 means stop at next NS step */
994 if ( (gs
.set
[eglsSTOPCOND
] < 0) ||
995 ( (gs
.set
[eglsSTOPCOND
] > 0) && (bNStList
|| ir
->nstlist
== 0) ) )
1000 /* Determine whether or not to update the Born radii if doing GB */
1001 bBornRadii
= bFirstStep
;
1002 if (ir
->implicit_solvent
&& (step
% ir
->nstgbradii
== 0))
1007 do_log
= do_per_step(step
, ir
->nstlog
) || bFirstStep
|| bLastStep
;
1008 do_verbose
= bVerbose
&&
1009 (step
% stepout
== 0 || bFirstStep
|| bLastStep
);
1011 if (bNS
&& !(bFirstStep
&& ir
->bContinuation
&& !bRerunMD
))
1015 bMasterState
= TRUE
;
1019 bMasterState
= FALSE
;
1020 /* Correct the new box if it is too skewed */
1021 if (DYNAMIC_BOX(*ir
))
1023 if (correct_box(fplog
, step
, state
->box
, graph
))
1025 bMasterState
= TRUE
;
1028 if (DOMAINDECOMP(cr
) && bMasterState
)
1030 dd_collect_state(cr
->dd
, state
, state_global
);
1034 if (DOMAINDECOMP(cr
))
1036 /* Repartition the domain decomposition */
1037 wallcycle_start(wcycle
, ewcDOMDEC
);
1038 dd_partition_system(fplog
, step
, cr
,
1039 bMasterState
, nstglobalcomm
,
1040 state_global
, top_global
, ir
,
1041 state
, &f
, mdatoms
, top
, fr
,
1042 vsite
, shellfc
, constr
,
1044 do_verbose
&& !bPMETuneRunning
);
1045 wallcycle_stop(wcycle
, ewcDOMDEC
);
1046 /* If using an iterative integrator, reallocate space to match the decomposition */
1050 if (MASTER(cr
) && do_log
&& !bFFscan
)
1052 print_ebin_header(fplog
, step
, t
, state
->lambda
[efptFEP
]); /* can we improve the information printed here? */
1055 if (ir
->efep
!= efepNO
)
1057 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
1060 if ((bRerunMD
&& rerun_fr
.bV
) || bExchanged
)
1063 /* We need the kinetic energy at minus the half step for determining
1064 * the full step kinetic energy and possibly for T-coupling.*/
1065 /* This may not be quite working correctly yet . . . . */
1066 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
1067 wcycle
, enerd
, NULL
, NULL
, NULL
, NULL
, mu_tot
,
1068 constr
, NULL
, FALSE
, state
->box
,
1069 top_global
, &pcurr
, top_global
->natoms
, &bSumEkinhOld
,
1070 CGLO_RERUNMD
| CGLO_GSTAT
| CGLO_TEMPERATURE
);
1072 clear_mat(force_vir
);
1074 /* Ionize the atoms if necessary */
1077 ionize(fplog
, oenv
, mdatoms
, top_global
, t
, ir
, state
->x
, state
->v
,
1078 mdatoms
->start
, mdatoms
->start
+mdatoms
->homenr
, state
->box
, cr
);
1081 /* Update force field in ffscan program */
1084 if (update_forcefield(fplog
,
1086 mdatoms
->nr
, state
->x
, state
->box
))
1094 GMX_MPE_LOG(ev_timestep2
);
1096 /* We write a checkpoint at this MD step when:
1097 * either at an NS step when we signalled through gs,
1098 * or at the last step (but not when we do not want confout),
1099 * but never at the first step or with rerun.
1101 bCPT
= (((gs
.set
[eglsCHKPT
] && (bNS
|| ir
->nstlist
== 0)) ||
1102 (bLastStep
&& (Flags
& MD_CONFOUT
))) &&
1103 step
> ir
->init_step
&& !bRerunMD
);
1106 gs
.set
[eglsCHKPT
] = 0;
1109 /* Determine the energy and pressure:
1110 * at nstcalcenergy steps and at energy output steps (set below).
1112 if (EI_VV(ir
->eI
) && (!bInitStep
))
1114 /* for vv, the first half of the integration actually corresponds
1115 to the previous step. bCalcEner is only required to be evaluated on the 'next' step,
1116 but the virial needs to be calculated on both the current step and the 'next' step. Future
1117 reorganization may be able to get rid of one of the bCalcVir=TRUE steps. */
1119 bCalcEner
= do_per_step(step
-1, ir
->nstcalcenergy
);
1120 bCalcVir
= bCalcEner
||
1121 (ir
->epc
!= epcNO
&& (do_per_step(step
, ir
->nstpcouple
) || do_per_step(step
-1, ir
->nstpcouple
)));
1125 bCalcEner
= do_per_step(step
, ir
->nstcalcenergy
);
1126 bCalcVir
= bCalcEner
||
1127 (ir
->epc
!= epcNO
&& do_per_step(step
, ir
->nstpcouple
));
1130 /* Do we need global communication ? */
1131 bGStat
= (bCalcVir
|| bCalcEner
|| bStopCM
||
1132 do_per_step(step
, nstglobalcomm
) || (bVV
&& IR_NVT_TROTTER(ir
) && do_per_step(step
-1, nstglobalcomm
)) ||
1133 (ir
->nstlist
== -1 && !bRerunMD
&& step
>= nlh
.step_nscheck
));
1135 do_ene
= (do_per_step(step
, ir
->nstenergy
) || bLastStep
);
1137 if (do_ene
|| do_log
)
1144 /* these CGLO_ options remain the same throughout the iteration */
1145 cglo_flags
= ((bRerunMD
? CGLO_RERUNMD
: 0) |
1146 (bGStat
? CGLO_GSTAT
: 0)
1149 force_flags
= (GMX_FORCE_STATECHANGED
|
1150 ((DYNAMIC_BOX(*ir
) || bRerunMD
) ? GMX_FORCE_DYNAMICBOX
: 0) |
1151 GMX_FORCE_ALLFORCES
|
1153 (bCalcVir
? GMX_FORCE_VIRIAL
: 0) |
1154 (bCalcEner
? GMX_FORCE_ENERGY
: 0) |
1155 (bDoFEP
? GMX_FORCE_DHDL
: 0)
1160 if (do_per_step(step
, ir
->nstcalclr
))
1162 force_flags
|= GMX_FORCE_DO_LR
;
1168 /* Now is the time to relax the shells */
1169 count
= relax_shell_flexcon(fplog
, cr
, bVerbose
, bFFscan
? step
+1 : step
,
1170 ir
, bNS
, force_flags
,
1171 bStopCM
, top
, top_global
,
1173 state
, f
, force_vir
, mdatoms
,
1174 nrnb
, wcycle
, graph
, groups
,
1175 shellfc
, fr
, bBornRadii
, t
, mu_tot
,
1176 state
->natoms
, &bConverged
, vsite
,
1187 /* The coordinates (x) are shifted (to get whole molecules)
1189 * This is parallellized as well, and does communication too.
1190 * Check comments in sim_util.c
1192 do_force(fplog
, cr
, ir
, step
, nrnb
, wcycle
, top
, top_global
, groups
,
1193 state
->box
, state
->x
, &state
->hist
,
1194 f
, force_vir
, mdatoms
, enerd
, fcd
,
1195 state
->lambda
, graph
,
1196 fr
, vsite
, mu_tot
, t
, outf
->fp_field
, ed
, bBornRadii
,
1197 (bNS
? GMX_FORCE_NS
: 0) | force_flags
);
1200 GMX_BARRIER(cr
->mpi_comm_mygroup
);
1204 mu_aver
= calc_mu_aver(cr
, state
->x
, mdatoms
->chargeA
,
1205 mu_tot
, &top_global
->mols
, mdatoms
, gnx
, grpindex
);
1208 if (bTCR
&& bFirstStep
)
1210 tcr
= init_coupling(fplog
, nfile
, fnm
, cr
, fr
, mdatoms
, &(top
->idef
));
1211 fprintf(fplog
, "Done init_coupling\n");
1215 if (bVV
&& !bStartingFromCpt
&& !bRerunMD
)
1216 /* ############### START FIRST UPDATE HALF-STEP FOR VV METHODS############### */
1218 if (ir
->eI
== eiVV
&& bInitStep
)
1220 /* if using velocity verlet with full time step Ekin,
1221 * take the first half step only to compute the
1222 * virial for the first step. From there,
1223 * revert back to the initial coordinates
1224 * so that the input is actually the initial step.
1226 copy_rvecn(state
->v
, cbuf
, 0, state
->natoms
); /* should make this better for parallelizing? */
1230 /* this is for NHC in the Ekin(t+dt/2) version of vv */
1231 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ1
);
1234 /* If we are using twin-range interactions where the long-range component
1235 * is only evaluated every nstcalclr>1 steps, we should do a special update
1236 * step to combine the long-range forces on these steps.
1237 * For nstcalclr=1 this is not done, since the forces would have been added
1238 * directly to the short-range forces already.
1240 * TODO Remove various aspects of VV+twin-range in master
1241 * branch, because VV integrators did not ever support
1242 * twin-range multiple time stepping with constraints.
1244 bUpdateDoLR
= (fr
->bTwinRange
&& do_per_step(step
, ir
->nstcalclr
));
1246 update_coords(fplog
, step
, ir
, mdatoms
, state
, fr
->bMolPBC
,
1247 f
, bUpdateDoLR
, fr
->f_twin
, bCalcVir
? &fr
->vir_twin_constr
: NULL
, fcd
,
1248 ekind
, M
, wcycle
, upd
, bInitStep
, etrtVELOCITY1
,
1249 cr
, nrnb
, constr
, &top
->idef
);
1251 if (bIterativeCase
&& do_per_step(step
-1, ir
->nstpcouple
) && !bInitStep
)
1253 gmx_iterate_init(&iterate
, TRUE
);
1255 /* for iterations, we save these vectors, as we will be self-consistently iterating
1258 /*#### UPDATE EXTENDED VARIABLES IN TROTTER FORMULATION */
1260 /* save the state */
1261 if (iterate
.bIterationActive
)
1263 copy_coupling_state(state
, bufstate
, ekind
, ekind_save
, &(ir
->opts
));
1266 bFirstIterate
= TRUE
;
1267 while (bFirstIterate
|| iterate
.bIterationActive
)
1269 if (iterate
.bIterationActive
)
1271 copy_coupling_state(bufstate
, state
, ekind_save
, ekind
, &(ir
->opts
));
1272 if (bFirstIterate
&& bTrotter
)
1274 /* The first time through, we need a decent first estimate
1275 of veta(t+dt) to compute the constraints. Do
1276 this by computing the box volume part of the
1277 trotter integration at this time. Nothing else
1278 should be changed by this routine here. If
1279 !(first time), we start with the previous value
1282 veta_save
= state
->veta
;
1283 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ0
);
1284 vetanew
= state
->veta
;
1285 state
->veta
= veta_save
;
1290 if (!bRerunMD
|| rerun_fr
.bV
|| bForceUpdate
) /* Why is rerun_fr.bV here? Unclear. */
1292 update_constraints(fplog
, step
, NULL
, ir
, ekind
, mdatoms
,
1293 state
, fr
->bMolPBC
, graph
, f
,
1294 &top
->idef
, shake_vir
, NULL
,
1295 cr
, nrnb
, wcycle
, upd
, constr
,
1296 bInitStep
, TRUE
, bCalcVir
, vetanew
);
1298 if (bCalcVir
&& bUpdateDoLR
&& ir
->nstcalclr
> 1)
1300 /* Correct the virial for multiple time stepping */
1301 m_sub(shake_vir
, fr
->vir_twin_constr
, shake_vir
);
1304 if (!bOK
&& !bFFscan
)
1306 gmx_fatal(FARGS
, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1312 /* Need to unshift here if a do_force has been
1313 called in the previous step */
1314 unshift_self(graph
, state
->box
, state
->x
);
1317 /* if VV, compute the pressure and constraints */
1318 /* For VV2, we strictly only need this if using pressure
1319 * control, but we really would like to have accurate pressures
1321 * Think about ways around this in the future?
1322 * For now, keep this choice in comments.
1324 /*bPres = (ir->eI==eiVV || IR_NPT_TROTTER(ir)); */
1325 /*bTemp = ((ir->eI==eiVV &&(!bInitStep)) || (ir->eI==eiVVAK && IR_NPT_TROTTER(ir)));*/
1327 bTemp
= ((ir
->eI
== eiVV
&& (!bInitStep
)) || (ir
->eI
== eiVVAK
));
1328 if (bCalcEner
&& ir
->eI
== eiVVAK
) /*MRS: 7/9/2010 -- this still doesn't fix it?*/
1330 bSumEkinhOld
= TRUE
;
1332 /* for vv, the first half of the integration actually corresponds to the previous step.
1333 So we need information from the last step in the first half of the integration */
1334 if (bGStat
|| do_per_step(step
-1, nstglobalcomm
))
1336 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
1337 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1338 constr
, NULL
, FALSE
, state
->box
,
1339 top_global
, &pcurr
, top_global
->natoms
, &bSumEkinhOld
,
1342 | (bTemp
? CGLO_TEMPERATURE
: 0)
1343 | (bPres
? CGLO_PRESSURE
: 0)
1344 | (bPres
? CGLO_CONSTRAINT
: 0)
1345 | ((iterate
.bIterationActive
) ? CGLO_ITERATE
: 0)
1346 | (bFirstIterate
? CGLO_FIRSTITERATE
: 0)
1349 /* explanation of above:
1350 a) We compute Ekin at the full time step
1351 if 1) we are using the AveVel Ekin, and it's not the
1352 initial step, or 2) if we are using AveEkin, but need the full
1353 time step kinetic energy for the pressure (always true now, since we want accurate statistics).
1354 b) If we are using EkinAveEkin for the kinetic energy for the temperature control, we still feed in
1355 EkinAveVel because it's needed for the pressure */
1357 /* temperature scaling and pressure scaling to produce the extended variables at t+dt */
1362 m_add(force_vir
, shake_vir
, total_vir
); /* we need the un-dispersion corrected total vir here */
1363 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ2
);
1370 /* We need the kinetic energy at minus the half step for determining
1371 * the full step kinetic energy and possibly for T-coupling.*/
1372 /* This may not be quite working correctly yet . . . . */
1373 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
1374 wcycle
, enerd
, NULL
, NULL
, NULL
, NULL
, mu_tot
,
1375 constr
, NULL
, FALSE
, state
->box
,
1376 top_global
, &pcurr
, top_global
->natoms
, &bSumEkinhOld
,
1377 CGLO_RERUNMD
| CGLO_GSTAT
| CGLO_TEMPERATURE
);
1382 if (iterate
.bIterationActive
&&
1383 done_iterating(cr
, fplog
, step
, &iterate
, bFirstIterate
,
1384 state
->veta
, &vetanew
))
1388 bFirstIterate
= FALSE
;
1391 if (bTrotter
&& !bInitStep
)
1393 copy_mat(shake_vir
, state
->svir_prev
);
1394 copy_mat(force_vir
, state
->fvir_prev
);
1395 if (IR_NVT_TROTTER(ir
) && ir
->eI
== eiVV
)
1397 /* update temperature and kinetic energy now that step is over - this is the v(t+dt) point */
1398 enerd
->term
[F_TEMP
] = sum_ekin(&(ir
->opts
), ekind
, NULL
, (ir
->eI
== eiVV
), FALSE
, FALSE
);
1399 enerd
->term
[F_EKIN
] = trace(ekind
->ekin
);
1402 /* if it's the initial step, we performed this first step just to get the constraint virial */
1403 if (bInitStep
&& ir
->eI
== eiVV
)
1405 copy_rvecn(cbuf
, state
->v
, 0, state
->natoms
);
1408 GMX_MPE_LOG(ev_timestep1
);
1411 /* MRS -- now done iterating -- compute the conserved quantity */
1414 saved_conserved_quantity
= compute_conserved_from_auxiliary(ir
, state
, &MassQ
);
1417 last_ekin
= enerd
->term
[F_EKIN
];
1419 if ((ir
->eDispCorr
!= edispcEnerPres
) && (ir
->eDispCorr
!= edispcAllEnerPres
))
1421 saved_conserved_quantity
-= enerd
->term
[F_DISPCORR
];
1423 /* sum up the foreign energy and dhdl terms for vv. currently done every step so that dhdl is correct in the .edr */
1426 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1430 /* ######## END FIRST UPDATE STEP ############## */
1431 /* ######## If doing VV, we now have v(dt) ###### */
1434 /* perform extended ensemble sampling in lambda - we don't
1435 actually move to the new state before outputting
1436 statistics, but if performing simulated tempering, we
1437 do update the velocities and the tau_t. */
1439 lamnew
= ExpandedEnsembleDynamics(fplog
, ir
, enerd
, state
, &MassQ
, state
->fep_state
, &state
->dfhist
, step
, mcrng
, state
->v
, mdatoms
);
1440 /* history is maintained in state->dfhist, but state_global is what is sent to trajectory and log output */
1441 copy_df_history(&state_global
->dfhist
,&state
->dfhist
);
1443 /* ################## START TRAJECTORY OUTPUT ################# */
1445 /* Now we have the energies and forces corresponding to the
1446 * coordinates at time t. We must output all of this before
1448 * for RerunMD t is read from input trajectory
1450 GMX_MPE_LOG(ev_output_start
);
1453 if (do_per_step(step
, ir
->nstxout
))
1455 mdof_flags
|= MDOF_X
;
1457 if (do_per_step(step
, ir
->nstvout
))
1459 mdof_flags
|= MDOF_V
;
1461 if (do_per_step(step
, ir
->nstfout
))
1463 mdof_flags
|= MDOF_F
;
1465 if (do_per_step(step
, ir
->nstxtcout
))
1467 mdof_flags
|= MDOF_XTC
;
1471 mdof_flags
|= MDOF_CPT
;
1475 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
1478 /* Enforce writing positions and velocities at end of run */
1479 mdof_flags
|= (MDOF_X
| MDOF_V
);
1485 fcReportProgress( ir
->nsteps
, step
);
1488 #if defined(__native_client__)
1489 fcCheckin(MASTER(cr
));
1492 /* sync bCPT and fc record-keeping */
1493 if (bCPT
&& MASTER(cr
))
1495 fcRequestCheckPoint();
1499 if (mdof_flags
!= 0)
1501 wallcycle_start(wcycle
, ewcTRAJ
);
1504 if (state
->flags
& (1<<estLD_RNG
))
1506 get_stochd_state(upd
, state
);
1508 if (state
->flags
& (1<<estMC_RNG
))
1510 get_mc_state(mcrng
, state
);
1516 state_global
->ekinstate
.bUpToDate
= FALSE
;
1520 update_ekinstate(&state_global
->ekinstate
, ekind
);
1521 state_global
->ekinstate
.bUpToDate
= TRUE
;
1523 update_energyhistory(&state_global
->enerhist
, mdebin
);
1526 write_traj(fplog
, cr
, outf
, mdof_flags
, top_global
,
1527 step
, t
, state
, state_global
, f
, f_global
, &n_xtc
, &x_xtc
);
1534 if (bLastStep
&& step_rel
== ir
->nsteps
&&
1535 (Flags
& MD_CONFOUT
) && MASTER(cr
) &&
1536 !bRerunMD
&& !bFFscan
)
1538 /* x and v have been collected in write_traj,
1539 * because a checkpoint file will always be written
1542 fprintf(stderr
, "\nWriting final coordinates.\n");
1545 /* Make molecules whole only for confout writing */
1546 do_pbc_mtop(fplog
, ir
->ePBC
, state
->box
, top_global
, state_global
->x
);
1548 write_sto_conf_mtop(ftp2fn(efSTO
, nfile
, fnm
),
1549 *top_global
->name
, top_global
,
1550 state_global
->x
, state_global
->v
,
1551 ir
->ePBC
, state
->box
);
1554 wallcycle_stop(wcycle
, ewcTRAJ
);
1556 GMX_MPE_LOG(ev_output_finish
);
1558 /* kludge -- virial is lost with restart for NPT control. Must restart */
1559 if (bStartingFromCpt
&& bVV
)
1561 copy_mat(state
->svir_prev
, shake_vir
);
1562 copy_mat(state
->fvir_prev
, force_vir
);
1564 /* ################## END TRAJECTORY OUTPUT ################ */
1566 /* Determine the wallclock run time up till now */
1567 run_time
= gmx_gettime() - (double)runtime
->real
;
1569 /* Check whether everything is still allright */
1570 if (((int)gmx_get_stop_condition() > handled_stop_condition
)
1571 #ifdef GMX_THREAD_MPI
1576 /* this is just make gs.sig compatible with the hack
1577 of sending signals around by MPI_Reduce with together with
1579 if (gmx_get_stop_condition() == gmx_stop_cond_next_ns
)
1581 gs
.sig
[eglsSTOPCOND
] = 1;
1583 if (gmx_get_stop_condition() == gmx_stop_cond_next
)
1585 gs
.sig
[eglsSTOPCOND
] = -1;
1587 /* < 0 means stop at next step, > 0 means stop at next NS step */
1591 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1592 gmx_get_signal_name(),
1593 gs
.sig
[eglsSTOPCOND
] == 1 ? "NS " : "");
1597 "\n\nReceived the %s signal, stopping at the next %sstep\n\n",
1598 gmx_get_signal_name(),
1599 gs
.sig
[eglsSTOPCOND
] == 1 ? "NS " : "");
1601 handled_stop_condition
= (int)gmx_get_stop_condition();
1603 else if (MASTER(cr
) && (bNS
|| ir
->nstlist
<= 0) &&
1604 (max_hours
> 0 && run_time
> max_hours
*60.0*60.0*0.99) &&
1605 gs
.sig
[eglsSTOPCOND
] == 0 && gs
.set
[eglsSTOPCOND
] == 0)
1607 /* Signal to terminate the run */
1608 gs
.sig
[eglsSTOPCOND
] = 1;
1611 fprintf(fplog
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step
, sbuf
), max_hours
*0.99);
1613 fprintf(stderr
, "\nStep %s: Run time exceeded %.3f hours, will terminate the run\n", gmx_step_str(step
, sbuf
), max_hours
*0.99);
1616 if (bResetCountersHalfMaxH
&& MASTER(cr
) &&
1617 run_time
> max_hours
*60.0*60.0*0.495)
1619 gs
.sig
[eglsRESETCOUNTERS
] = 1;
1622 if (ir
->nstlist
== -1 && !bRerunMD
)
1624 /* When bGStatEveryStep=FALSE, global_stat is only called
1625 * when we check the atom displacements, not at NS steps.
1626 * This means that also the bonded interaction count check is not
1627 * performed immediately after NS. Therefore a few MD steps could
1628 * be performed with missing interactions.
1629 * But wrong energies are never written to file,
1630 * since energies are only written after global_stat
1633 if (step
>= nlh
.step_nscheck
)
1635 nlh
.nabnsb
= natoms_beyond_ns_buffer(ir
, fr
, &top
->cgs
,
1636 nlh
.scale_tot
, state
->x
);
1640 /* This is not necessarily true,
1641 * but step_nscheck is determined quite conservatively.
1647 /* In parallel we only have to check for checkpointing in steps
1648 * where we do global communication,
1649 * otherwise the other nodes don't know.
1651 if (MASTER(cr
) && ((bGStat
|| !PAR(cr
)) &&
1654 run_time
>= nchkpt
*cpt_period
*60.0)) &&
1655 gs
.set
[eglsCHKPT
] == 0)
1657 gs
.sig
[eglsCHKPT
] = 1;
1660 /* at the start of step, randomize or scale the velocities (trotter done elsewhere) */
1665 update_tcouple(fplog
, step
, ir
, state
, ekind
, wcycle
, upd
, &MassQ
, mdatoms
);
1667 if (ETC_ANDERSEN(ir
->etc
)) /* keep this outside of update_tcouple because of the extra info required to pass */
1669 gmx_bool bIfRandomize
;
1670 bIfRandomize
= update_randomize_velocities(ir
, step
, mdatoms
, state
, upd
, &top
->idef
, constr
, DOMAINDECOMP(cr
));
1671 /* if we have constraints, we have to remove the kinetic energy parallel to the bonds */
1672 if (constr
&& bIfRandomize
)
1674 update_constraints(fplog
, step
, NULL
, ir
, ekind
, mdatoms
,
1675 state
, fr
->bMolPBC
, graph
, f
,
1676 &top
->idef
, tmp_vir
, NULL
,
1677 cr
, nrnb
, wcycle
, upd
, constr
,
1678 bInitStep
, TRUE
, bCalcVir
, vetanew
);
1683 if (bIterativeCase
&& do_per_step(step
, ir
->nstpcouple
))
1685 gmx_iterate_init(&iterate
, TRUE
);
1686 /* for iterations, we save these vectors, as we will be redoing the calculations */
1687 copy_coupling_state(state
, bufstate
, ekind
, ekind_save
, &(ir
->opts
));
1690 bFirstIterate
= TRUE
;
1691 while (bFirstIterate
|| iterate
.bIterationActive
)
1693 /* We now restore these vectors to redo the calculation with improved extended variables */
1694 if (iterate
.bIterationActive
)
1696 copy_coupling_state(bufstate
, state
, ekind_save
, ekind
, &(ir
->opts
));
1699 /* We make the decision to break or not -after- the calculation of Ekin and Pressure,
1700 so scroll down for that logic */
1702 /* ######### START SECOND UPDATE STEP ################# */
1703 GMX_MPE_LOG(ev_update_start
);
1704 /* Box is changed in update() when we do pressure coupling,
1705 * but we should still use the old box for energy corrections and when
1706 * writing it to the energy file, so it matches the trajectory files for
1707 * the same timestep above. Make a copy in a separate array.
1709 copy_mat(state
->box
, lastbox
);
1714 if (!(bRerunMD
&& !rerun_fr
.bV
&& !bForceUpdate
))
1716 wallcycle_start(wcycle
, ewcUPDATE
);
1717 /* UPDATE PRESSURE VARIABLES IN TROTTER FORMULATION WITH CONSTRAINTS */
1720 if (iterate
.bIterationActive
)
1728 /* we use a new value of scalevir to converge the iterations faster */
1729 scalevir
= tracevir
/trace(shake_vir
);
1731 msmul(shake_vir
, scalevir
, shake_vir
);
1732 m_add(force_vir
, shake_vir
, total_vir
);
1733 clear_mat(shake_vir
);
1735 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ3
);
1736 /* We can only do Berendsen coupling after we have summed
1737 * the kinetic energy or virial. Since the happens
1738 * in global_state after update, we should only do it at
1739 * step % nstlist = 1 with bGStatEveryStep=FALSE.
1744 update_tcouple(fplog
, step
, ir
, state
, ekind
, wcycle
, upd
, &MassQ
, mdatoms
);
1745 update_pcouple(fplog
, step
, ir
, state
, pcoupl_mu
, M
, wcycle
,
1751 bUpdateDoLR
= (fr
->bTwinRange
&& do_per_step(step
, ir
->nstcalclr
));
1753 /* velocity half-step update */
1754 update_coords(fplog
, step
, ir
, mdatoms
, state
, fr
->bMolPBC
, f
,
1755 bUpdateDoLR
, fr
->f_twin
, bCalcVir
? &fr
->vir_twin_constr
: NULL
, fcd
,
1756 ekind
, M
, wcycle
, upd
, FALSE
, etrtVELOCITY2
,
1757 cr
, nrnb
, constr
, &top
->idef
);
1760 /* Above, initialize just copies ekinh into ekin,
1761 * it doesn't copy position (for VV),
1762 * and entire integrator for MD.
1765 if (ir
->eI
== eiVVAK
)
1767 copy_rvecn(state
->x
, cbuf
, 0, state
->natoms
);
1769 bUpdateDoLR
= (fr
->bTwinRange
&& do_per_step(step
, ir
->nstcalclr
));
1771 update_coords(fplog
, step
, ir
, mdatoms
, state
, fr
->bMolPBC
, f
,
1772 bUpdateDoLR
, fr
->f_twin
, bCalcVir
? &fr
->vir_twin_constr
: NULL
, fcd
,
1773 ekind
, M
, wcycle
, upd
, bInitStep
, etrtPOSITION
, cr
, nrnb
, constr
, &top
->idef
);
1774 wallcycle_stop(wcycle
, ewcUPDATE
);
1776 update_constraints(fplog
, step
, &dvdl_constr
, ir
, ekind
, mdatoms
, state
,
1777 fr
->bMolPBC
, graph
, f
,
1778 &top
->idef
, shake_vir
, force_vir
,
1779 cr
, nrnb
, wcycle
, upd
, constr
,
1780 bInitStep
, FALSE
, bCalcVir
, state
->veta
);
1782 if (bCalcVir
&& bUpdateDoLR
&& ir
->nstcalclr
> 1)
1784 /* Correct the virial for multiple time stepping */
1785 m_sub(shake_vir
, fr
->vir_twin_constr
, shake_vir
);
1788 if (ir
->eI
== eiVVAK
)
1790 /* erase F_EKIN and F_TEMP here? */
1791 /* just compute the kinetic energy at the half step to perform a trotter step */
1792 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
1793 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1794 constr
, NULL
, FALSE
, lastbox
,
1795 top_global
, &pcurr
, top_global
->natoms
, &bSumEkinhOld
,
1796 cglo_flags
| CGLO_TEMPERATURE
1798 wallcycle_start(wcycle
, ewcUPDATE
);
1799 trotter_update(ir
, step
, ekind
, enerd
, state
, total_vir
, mdatoms
, &MassQ
, trotter_seq
, ettTSEQ4
);
1800 /* now we know the scaling, we can compute the positions again again */
1801 copy_rvecn(cbuf
, state
->x
, 0, state
->natoms
);
1803 bUpdateDoLR
= (fr
->bTwinRange
&& do_per_step(step
, ir
->nstcalclr
));
1805 update_coords(fplog
, step
, ir
, mdatoms
, state
, fr
->bMolPBC
, f
,
1806 bUpdateDoLR
, fr
->f_twin
, bCalcVir
? &fr
->vir_twin_constr
: NULL
, fcd
,
1807 ekind
, M
, wcycle
, upd
, bInitStep
, etrtPOSITION
, cr
, nrnb
, constr
, &top
->idef
);
1808 wallcycle_stop(wcycle
, ewcUPDATE
);
1810 /* do we need an extra constraint here? just need to copy out of state->v to upd->xp? */
1811 /* are the small terms in the shake_vir here due
1812 * to numerical errors, or are they important
1813 * physically? I'm thinking they are just errors, but not completely sure.
1814 * For now, will call without actually constraining, constr=NULL*/
1815 update_constraints(fplog
, step
, NULL
, ir
, ekind
, mdatoms
,
1816 state
, fr
->bMolPBC
, graph
, f
,
1817 &top
->idef
, tmp_vir
, force_vir
,
1818 cr
, nrnb
, wcycle
, upd
, NULL
,
1819 bInitStep
, FALSE
, bCalcVir
,
1822 if (!bOK
&& !bFFscan
)
1824 gmx_fatal(FARGS
, "Constraint error: Shake, Lincs or Settle could not solve the constrains");
1827 if (fr
->bSepDVDL
&& fplog
&& do_log
)
1829 fprintf(fplog
, sepdvdlformat
, "Constraint dV/dl", 0.0, dvdl_constr
);
1833 /* this factor or 2 correction is necessary
1834 because half of the constraint force is removed
1835 in the vv step, so we have to double it. See
1836 the Redmine issue #1255. It is not yet clear
1837 if the factor of 2 is exact, or just a very
1838 good approximation, and this will be
1839 investigated. The next step is to see if this
1840 can be done adding a dhdl contribution from the
1841 rattle step, but this is somewhat more
1842 complicated with the current code. Will be
1843 investigated, hopefully for 4.6.3. However,
1844 this current solution is much better than
1845 having it completely wrong.
1847 enerd
->term
[F_DVDL_CONSTR
] += 2*dvdl_constr
;
1851 enerd
->term
[F_DVDL_CONSTR
] += dvdl_constr
;
1856 /* Need to unshift here */
1857 unshift_self(graph
, state
->box
, state
->x
);
1860 GMX_BARRIER(cr
->mpi_comm_mygroup
);
1861 GMX_MPE_LOG(ev_update_finish
);
1865 wallcycle_start(wcycle
, ewcVSITECONSTR
);
1868 shift_self(graph
, state
->box
, state
->x
);
1870 construct_vsites(fplog
, vsite
, state
->x
, nrnb
, ir
->delta_t
, state
->v
,
1871 top
->idef
.iparams
, top
->idef
.il
,
1872 fr
->ePBC
, fr
->bMolPBC
, graph
, cr
, state
->box
);
1876 unshift_self(graph
, state
->box
, state
->x
);
1878 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
1881 /* ############## IF NOT VV, Calculate globals HERE, also iterate constraints ############ */
1882 /* With Leap-Frog we can skip compute_globals at
1883 * non-communication steps, but we need to calculate
1884 * the kinetic energy one step before communication.
1886 if (bGStat
|| (!EI_VV(ir
->eI
) && do_per_step(step
+1, nstglobalcomm
)))
1888 if (ir
->nstlist
== -1 && bFirstIterate
)
1890 gs
.sig
[eglsNABNSB
] = nlh
.nabnsb
;
1892 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, state_global
, mdatoms
, nrnb
, vcm
,
1893 wcycle
, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
1895 bFirstIterate
? &gs
: NULL
,
1896 (step_rel
% gs
.nstms
== 0) &&
1897 (multisim_nsteps
< 0 || (step_rel
< multisim_nsteps
)),
1899 top_global
, &pcurr
, top_global
->natoms
, &bSumEkinhOld
,
1901 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_ENERGY
: 0)
1902 | (!EI_VV(ir
->eI
) && bStopCM
? CGLO_STOPCM
: 0)
1903 | (!EI_VV(ir
->eI
) ? CGLO_TEMPERATURE
: 0)
1904 | (!EI_VV(ir
->eI
) || bRerunMD
? CGLO_PRESSURE
: 0)
1905 | (iterate
.bIterationActive
? CGLO_ITERATE
: 0)
1906 | (bFirstIterate
? CGLO_FIRSTITERATE
: 0)
1909 if (ir
->nstlist
== -1 && bFirstIterate
)
1911 nlh
.nabnsb
= gs
.set
[eglsNABNSB
];
1912 gs
.set
[eglsNABNSB
] = 0;
1915 /* bIterate is set to keep it from eliminating the old ekin kinetic energy terms */
1916 /* ############# END CALC EKIN AND PRESSURE ################# */
1918 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
1919 the virial that should probably be addressed eventually. state->veta has better properies,
1920 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
1921 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
1923 if (iterate
.bIterationActive
&&
1924 done_iterating(cr
, fplog
, step
, &iterate
, bFirstIterate
,
1925 trace(shake_vir
), &tracevir
))
1929 bFirstIterate
= FALSE
;
1932 if (!bVV
|| bRerunMD
)
1934 /* sum up the foreign energy and dhdl terms for md and sd. currently done every step so that dhdl is correct in the .edr */
1935 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
1937 update_box(fplog
, step
, ir
, mdatoms
, state
, graph
, f
,
1938 ir
->nstlist
== -1 ? &nlh
.scale_tot
: NULL
, pcoupl_mu
, nrnb
, wcycle
, upd
, bInitStep
, FALSE
);
1940 /* ################# END UPDATE STEP 2 ################# */
1941 /* #### We now have r(t+dt) and v(t+dt/2) ############# */
1943 /* The coordinates (x) were unshifted in update */
1944 if (bFFscan
&& (shellfc
== NULL
|| bConverged
))
1946 if (print_forcefield(fplog
, enerd
->term
, mdatoms
->homenr
,
1948 &(top_global
->mols
), mdatoms
->massT
, pres
))
1952 fprintf(stderr
, "\n");
1958 /* We will not sum ekinh_old,
1959 * so signal that we still have to do it.
1961 bSumEkinhOld
= TRUE
;
1966 /* Only do GCT when the relaxation of shells (minimization) has converged,
1967 * otherwise we might be coupling to bogus energies.
1968 * In parallel we must always do this, because the other sims might
1972 /* Since this is called with the new coordinates state->x, I assume
1973 * we want the new box state->box too. / EL 20040121
1975 do_coupling(fplog
, oenv
, nfile
, fnm
, tcr
, t
, step
, enerd
->term
, fr
,
1977 mdatoms
, &(top
->idef
), mu_aver
,
1978 top_global
->mols
.nr
, cr
,
1979 state
->box
, total_vir
, pres
,
1980 mu_tot
, state
->x
, f
, bConverged
);
1984 /* ######### BEGIN PREPARING EDR OUTPUT ########### */
1986 /* use the directly determined last velocity, not actually the averaged half steps */
1987 if (bTrotter
&& ir
->eI
== eiVV
)
1989 enerd
->term
[F_EKIN
] = last_ekin
;
1991 enerd
->term
[F_ETOT
] = enerd
->term
[F_EPOT
] + enerd
->term
[F_EKIN
];
1995 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + saved_conserved_quantity
;
1999 enerd
->term
[F_ECONSERVED
] = enerd
->term
[F_ETOT
] + compute_conserved_from_auxiliary(ir
, state
, &MassQ
);
2001 /* Check for excessively large energies */
2005 real etot_max
= 1e200
;
2007 real etot_max
= 1e30
;
2009 if (fabs(enerd
->term
[F_ETOT
]) > etot_max
)
2011 fprintf(stderr
, "Energy too large (%g), giving up\n",
2012 enerd
->term
[F_ETOT
]);
2015 /* ######### END PREPARING EDR OUTPUT ########### */
2017 /* Time for performance */
2018 if (((step
% stepout
) == 0) || bLastStep
)
2020 runtime_upd_proc(runtime
);
2026 gmx_bool do_dr
, do_or
;
2028 if (fplog
&& do_log
&& bDoExpanded
)
2030 /* only needed if doing expanded ensemble */
2031 PrintFreeEnergyInfoToFile(fplog
, ir
->fepvals
, ir
->expandedvals
, ir
->bSimTemp
? ir
->simtempvals
: NULL
,
2032 &state_global
->dfhist
, state
->fep_state
, ir
->nstlog
, step
);
2034 if (!(bStartingFromCpt
&& (EI_VV(ir
->eI
))))
2038 upd_mdebin(mdebin
, bDoDHDL
, TRUE
,
2039 t
, mdatoms
->tmass
, enerd
, state
,
2040 ir
->fepvals
, ir
->expandedvals
, lastbox
,
2041 shake_vir
, force_vir
, total_vir
, pres
,
2042 ekind
, mu_tot
, constr
);
2046 upd_mdebin_step(mdebin
);
2049 do_dr
= do_per_step(step
, ir
->nstdisreout
);
2050 do_or
= do_per_step(step
, ir
->nstorireout
);
2052 print_ebin(outf
->fp_ene
, do_ene
, do_dr
, do_or
, do_log
? fplog
: NULL
,
2054 eprNORMAL
, bCompact
, mdebin
, fcd
, groups
, &(ir
->opts
));
2056 if (ir
->ePull
!= epullNO
)
2058 pull_print_output(ir
->pull
, step
, t
);
2061 if (do_per_step(step
, ir
->nstlog
))
2063 if (fflush(fplog
) != 0)
2065 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
2071 /* Have to do this part _after_ outputting the logfile and the edr file */
2072 /* Gets written into the state at the beginning of next loop*/
2073 state
->fep_state
= lamnew
;
2076 /* Remaining runtime */
2077 if (MULTIMASTER(cr
) && (do_verbose
|| gmx_got_usr_signal()) && !bPMETuneRunning
)
2081 fprintf(stderr
, "\n");
2083 print_time(stderr
, runtime
, step
, ir
, cr
);
2086 /* Replica exchange */
2088 if ((repl_ex_nst
> 0) && (step
> 0) && !bLastStep
&&
2089 do_per_step(step
, repl_ex_nst
))
2091 bExchanged
= replica_exchange(fplog
, cr
, repl_ex
,
2092 state_global
, enerd
,
2095 if (bExchanged
&& DOMAINDECOMP(cr
))
2097 dd_partition_system(fplog
, step
, cr
, TRUE
, 1,
2098 state_global
, top_global
, ir
,
2099 state
, &f
, mdatoms
, top
, fr
,
2100 vsite
, shellfc
, constr
,
2101 nrnb
, wcycle
, FALSE
);
2107 bStartingFromCpt
= FALSE
;
2109 /* ####### SET VARIABLES FOR NEXT ITERATION IF THEY STILL NEED IT ###### */
2110 /* With all integrators, except VV, we need to retain the pressure
2111 * at the current step for coupling at the next step.
2113 if ((state
->flags
& (1<<estPRES_PREV
)) &&
2115 (ir
->nstpcouple
> 0 && step
% ir
->nstpcouple
== 0)))
2117 /* Store the pressure in t_state for pressure coupling
2118 * at the next MD step.
2120 copy_mat(pres
, state
->pres_prev
);
2123 /* ####### END SET VARIABLES FOR NEXT ITERATION ###### */
2125 if ( (membed
!= NULL
) && (!bLastStep
) )
2127 rescale_membed(step_rel
, membed
, state_global
->x
);
2134 /* read next frame from input trajectory */
2135 bNotLastFrame
= read_next_frame(oenv
, status
, &rerun_fr
);
2140 rerun_parallel_comm(cr
, &rerun_fr
, &bNotLastFrame
);
2144 if (!bRerunMD
|| !rerun_fr
.bStep
)
2146 /* increase the MD step number */
2151 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
2152 if (DOMAINDECOMP(cr
) && wcycle
)
2154 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
2157 if (bPMETuneRunning
|| bPMETuneTry
)
2159 /* PME grid + cut-off optimization with GPUs or PME nodes */
2161 /* Count the total cycles over the last steps */
2162 cycles_pmes
+= cycles
;
2164 /* We can only switch cut-off at NS steps */
2165 if (step
% ir
->nstlist
== 0)
2167 /* PME grid + cut-off optimization with GPUs or PME nodes */
2170 if (DDMASTER(cr
->dd
))
2172 /* PME node load is too high, start tuning */
2173 bPMETuneRunning
= (dd_pme_f_ratio(cr
->dd
) >= 1.05);
2175 dd_bcast(cr
->dd
, sizeof(gmx_bool
), &bPMETuneRunning
);
2177 if (bPMETuneRunning
|| step_rel
> ir
->nstlist
*50)
2179 bPMETuneTry
= FALSE
;
2182 if (bPMETuneRunning
)
2184 /* init_step might not be a multiple of nstlist,
2185 * but the first cycle is always skipped anyhow.
2188 pme_load_balance(pme_loadbal
, cr
,
2189 (bVerbose
&& MASTER(cr
)) ? stderr
: NULL
,
2191 ir
, state
, cycles_pmes
,
2192 fr
->ic
, fr
->nbv
, &fr
->pmedata
,
2195 /* Update constants in forcerec/inputrec to keep them in sync with fr->ic */
2196 fr
->ewaldcoeff
= fr
->ic
->ewaldcoeff
;
2197 fr
->rlist
= fr
->ic
->rlist
;
2198 fr
->rlistlong
= fr
->ic
->rlistlong
;
2199 fr
->rcoulomb
= fr
->ic
->rcoulomb
;
2200 fr
->rvdw
= fr
->ic
->rvdw
;
2206 if (step_rel
== wcycle_get_reset_counters(wcycle
) ||
2207 gs
.set
[eglsRESETCOUNTERS
] != 0)
2209 /* Reset all the counters related to performance over the run */
2210 reset_all_counters(fplog
, cr
, step
, &step_rel
, ir
, wcycle
, nrnb
, runtime
,
2211 fr
->nbv
!= NULL
&& fr
->nbv
->bUseGPU
? fr
->nbv
->cu_nbv
: NULL
);
2212 wcycle_set_reset_counters(wcycle
, -1);
2213 if (!(cr
->duty
& DUTY_PME
))
2215 /* Tell our PME node to reset its counters */
2216 gmx_pme_send_resetcounters(cr
, step
);
2218 /* Correct max_hours for the elapsed time */
2219 max_hours
-= run_time
/(60.0*60.0);
2220 bResetCountersHalfMaxH
= FALSE
;
2221 gs
.set
[eglsRESETCOUNTERS
] = 0;
2225 /* End of main MD loop */
2229 runtime_end(runtime
);
2231 if (bRerunMD
&& MASTER(cr
))
2236 if (!(cr
->duty
& DUTY_PME
))
2238 /* Tell the PME only node to finish */
2239 gmx_pme_send_finish(cr
);
2244 if (ir
->nstcalcenergy
> 0 && !bRerunMD
)
2246 print_ebin(outf
->fp_ene
, FALSE
, FALSE
, FALSE
, fplog
, step
, t
,
2247 eprAVER
, FALSE
, mdebin
, fcd
, groups
, &(ir
->opts
));
2255 if (ir
->nstlist
== -1 && nlh
.nns
> 0 && fplog
)
2257 fprintf(fplog
, "Average neighborlist lifetime: %.1f steps, std.dev.: %.1f steps\n", nlh
.s1
/nlh
.nns
, sqrt(nlh
.s2
/nlh
.nns
- sqr(nlh
.s1
/nlh
.nns
)));
2258 fprintf(fplog
, "Average number of atoms that crossed the half buffer length: %.1f\n\n", nlh
.ab
/nlh
.nns
);
2261 if (pme_loadbal
!= NULL
)
2263 pme_loadbal_done(pme_loadbal
, cr
, fplog
,
2264 fr
->nbv
!= NULL
&& fr
->nbv
->bUseGPU
);
2267 if (shellfc
&& fplog
)
2269 fprintf(fplog
, "Fraction of iterations that converged: %.2f %%\n",
2270 (nconverged
*100.0)/step_rel
);
2271 fprintf(fplog
, "Average number of force evaluations per MD step: %.2f\n\n",
2275 if (repl_ex_nst
> 0 && MASTER(cr
))
2277 print_replica_exchange_statistics(fplog
, repl_ex
);
2280 runtime
->nsteps_done
= step_rel
;