From: Mark Abraham Date: Tue, 8 Jul 2014 20:39:52 +0000 (+0200) Subject: Fix gmx tune_pme with LJ-PME X-Git-Url: http://biod.pnpi.spb.ru/gitweb/?a=commitdiff_plain;h=660199b14adac6a28adb0a3a7d21d1f88e62c99e;p=alexxy%2Fgromacs.git Fix gmx tune_pme with LJ-PME Extend the range of vdwtype for which tune_pme will scale rvdw. This makes tune_pme worth using with LJ-PME. Change-Id: Iec702ae984cd062d47970380c0ca41f82e4c31d2 --- diff --git a/src/gromacs/gmxana/gmx_tune_pme.c b/src/gromacs/gmxana/gmx_tune_pme.c index 361c6a82b3..112e872732 100644 --- a/src/gromacs/gmxana/gmx_tune_pme.c +++ b/src/gromacs/gmxana/gmx_tune_pme.c @@ -809,6 +809,11 @@ static void modify_PMEsettings( sfree(ir); } +static gmx_bool can_scale_rvdw(int vdwtype) +{ + return (evdwCUT == vdwtype || + evdwPME == vdwtype); +} #define EPME_SWITCHED(e) ((e) == eelPMESWITCH || (e) == eelPMEUSERSWITCH) @@ -958,7 +963,7 @@ static void make_benchmark_tprs( fprintf(fp, " No. scaling rcoulomb"); fprintf(fp, " nkx nky nkz"); fprintf(fp, " spacing"); - if (evdwCUT == ir->vdwtype) + if (can_scale_rvdw(ir->vdwtype)) { fprintf(fp, " rvdw"); } @@ -1015,11 +1020,14 @@ static void make_benchmark_tprs( ir->rlist = ir->rcoulomb + nlist_buffer; } - if (bScaleRvdw && evdwCUT == ir->vdwtype) + if (bScaleRvdw && can_scale_rvdw(ir->vdwtype)) { - if (ecutsVERLET == ir->cutoff_scheme) + if (ecutsVERLET == ir->cutoff_scheme || + evdwPME == ir->vdwtype) { - /* With Verlet, the van der Waals radius must always equal the Coulomb radius */ + /* With either the Verlet cutoff-scheme or LJ-PME, + the van der Waals radius must always equal the + Coulomb radius */ ir->rvdw = ir->rcoulomb; } else @@ -1067,7 +1075,7 @@ static void make_benchmark_tprs( fprintf(fp, "%4d%10f%10f", j, fac, ir->rcoulomb); fprintf(fp, "%5d%5d%5d", ir->nkx, ir->nky, ir->nkz); fprintf(fp, " %9f ", info->fsx[j]); - if (evdwCUT == ir->vdwtype) + if (can_scale_rvdw(ir->vdwtype)) { fprintf(fp, "%10f", ir->rvdw); }