From ec74f5259425aa6df35fbca95c018124d156f14a Mon Sep 17 00:00:00 2001 From: Berk Hess Date: Tue, 10 Jun 2014 14:38:04 +0200 Subject: [PATCH] Decoupled repl_ex_nst from nstcalcenergy The replica exchange frequency is automatically changed by mdrun to a multiple of nstcalcenergy, which is annoying. It turns out that it doesn't need to be a multiple, so this changing has been removed. Fixes #1494 Change-Id: I00833f92fd468924f61879aff8b7c85fe79d3c2e --- src/programs/mdrun/md.c | 20 +++++++------------- 1 file changed, 7 insertions(+), 13 deletions(-) diff --git a/src/programs/mdrun/md.c b/src/programs/mdrun/md.c index 60b01edd93..b7dab3cf20 100644 --- a/src/programs/mdrun/md.c +++ b/src/programs/mdrun/md.c @@ -191,7 +191,7 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[], int count, nconverged = 0; real timestep = 0; double tcount = 0; - gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bExchanged, bNeedRepartition; + gmx_bool bConverged = TRUE, bOK, bSumEkinhOld, bDoReplEx, bExchanged, bNeedRepartition; gmx_bool bAppend; gmx_bool bResetCountersHalfMaxH = FALSE; gmx_bool bVV, bIterativeCase, bFirstIterate, bTemp, bPres, bTrotter; @@ -465,15 +465,6 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[], set_constraints(constr, top, ir, mdatoms, cr); } - if (repl_ex_nst > 0) - { - /* We need to be sure replica exchange can only occur - * when the energies are current */ - check_nst_param(fplog, cr, "nstcalcenergy", ir->nstcalcenergy, - "repl_ex_nst", &repl_ex_nst); - /* This check needs to happen before inter-simulation - * signals are initialized, too */ - } if (repl_ex_nst > 0 && MASTER(cr)) { repl_ex = init_replica_exchange(fplog, cr->ms, state_global, ir, @@ -730,6 +721,7 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[], bStartingFromCpt = (Flags & MD_STARTFROMCPT) && bInitStep; bLastStep = FALSE; bSumEkinhOld = FALSE; + bDoReplEx = FALSE; bExchanged = FALSE; bNeedRepartition = FALSE; @@ -792,6 +784,9 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[], && (ir->bExpanded) && (step > 0) && (!bStartingFromCpt)); } + bDoReplEx = ((repl_ex_nst > 0) && (step > 0) && !bLastStep && + do_per_step(step, repl_ex_nst)); + if (bSimAnn) { update_annealing_target_temp(&(ir->opts), t); @@ -1021,7 +1016,7 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[], do_ene = (do_per_step(step, ir->nstenergy) || bLastStep); - if (do_ene || do_log) + if (do_ene || do_log || bDoReplEx) { bCalcVir = TRUE; bCalcEner = TRUE; @@ -1812,8 +1807,7 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[], /* Replica exchange */ bExchanged = FALSE; - if ((repl_ex_nst > 0) && (step > 0) && !bLastStep && - do_per_step(step, repl_ex_nst)) + if (bDoReplEx) { bExchanged = replica_exchange(fplog, cr, repl_ex, state_global, enerd, -- 2.22.0