From: Magnus Lundborg Date: Tue, 25 Sep 2018 08:38:03 +0000 (+0200) Subject: Add margin as input to pullGroupObeysPbcRestrictions X-Git-Url: http://biod.pnpi.spb.ru/gitweb/?a=commitdiff_plain;h=0d5de9982640d82142dd6b9e18bc5656e7357008;p=alexxy%2Fgromacs.git Add margin as input to pullGroupObeysPbcRestrictions Make it possible to specify the margin when checking that a pull group is not too large (compared to the PBC dimensions). Change-Id: I1a856ae6078c085693ed7add4cf765b5a49abc2d --- diff --git a/src/gromacs/gmxpreprocess/readpull.cpp b/src/gromacs/gmxpreprocess/readpull.cpp index a6de96cb82..db4f6a601c 100644 --- a/src/gromacs/gmxpreprocess/readpull.cpp +++ b/src/gromacs/gmxpreprocess/readpull.cpp @@ -522,7 +522,7 @@ pull_t *set_pull_init(t_inputrec *ir, const gmx_mtop_t *mtop, pull_calc_coms(nullptr, pull_work, md, &pbc, t_start, x, nullptr); - int groupThatFailsPbc = pullCheckPbcWithinGroups(*pull_work, x, pbc); + int groupThatFailsPbc = pullCheckPbcWithinGroups(*pull_work, x, pbc, c_pullGroupPbcMargin); if (groupThatFailsPbc >= 0) { char buf[STRLEN]; diff --git a/src/gromacs/pulling/pull.h b/src/gromacs/pulling/pull.h index ccb7cb8ab0..9ea3eefe20 100644 --- a/src/gromacs/pulling/pull.h +++ b/src/gromacs/pulling/pull.h @@ -288,7 +288,7 @@ void pull_calc_coms(const t_commrec *cr, const rvec x[], rvec *xp); -/*! \brief Margin for checking PBC distances compared to half the box size in pullCheckPbcWithinGroups() */ +/*! \brief Margin for checking pull group PBC distances compared to half the box size */ static constexpr real c_pullGroupPbcMargin = 0.9; /*! \brief Checks whether all groups that use a reference atom are within PBC restrictions @@ -297,19 +297,22 @@ static constexpr real c_pullGroupPbcMargin = 0.9; * atoms within half the box size from the PBC atom. The box size is used * per dimension for rectangular boxes, but can be a combination of * dimensions for triclinic boxes, depending on which dimensions are - * involved in the pull coordinates a group is involved in. + * involved in the pull coordinates a group is involved in. A margin is specified + * to ensure that atoms are not too close to the maximum distance. * * Should be called without MPI parallelization and after pull_calc_coms() * has been called at least once. * - * \param[in] pull The pull data structure - * \param[in] x The coordinates - * \param[in] pbc Information struct about periodicity + * \param[in] pull The pull data structure + * \param[in] x The coordinates + * \param[in] pbc Information struct about periodicity + * \param[in] pbcMargin The minimum margin (as a fraction) to half the box size * \returns -1 when all groups obey PBC or the first group index that fails PBC */ int pullCheckPbcWithinGroups(const pull_t &pull, const rvec *x, - const t_pbc &pbc); + const t_pbc &pbc, + real pbcMargin); /*! \brief Returns if we have pull coordinates with potential pulling. * diff --git a/src/gromacs/pulling/pullutil.cpp b/src/gromacs/pulling/pullutil.cpp index 5395ca0595..8280f6d61f 100644 --- a/src/gromacs/pulling/pullutil.cpp +++ b/src/gromacs/pulling/pullutil.cpp @@ -808,7 +808,8 @@ static bool pullGroupObeysPbcRestrictions(const pull_group_work_t &group, const BoolVec &dimUsed, const rvec *x, const t_pbc &pbc, - const gmx::RVec &x_pbc) + const gmx::RVec &x_pbc, + const real pbcMargin) { /* Determine which dimensions are relevant for PBC */ BoolVec dimUsesPbc = { false, false, false }; @@ -838,7 +839,7 @@ static bool pullGroupObeysPbcRestrictions(const pull_group_work_t &group, /* Use margins for dimensions independently */ for (int d = 0; d < pbc.ndim_ePBC; d++) { - marginPerDim[d] = c_pullGroupPbcMargin*pbc.hbox_diag[d]; + marginPerDim[d] = pbcMargin*pbc.hbox_diag[d]; } } else @@ -848,7 +849,7 @@ static bool pullGroupObeysPbcRestrictions(const pull_group_work_t &group, { if (dimUsesPbc[d]) { - marginDistance2 += c_pullGroupPbcMargin*gmx::square(0.5)*norm2(pbc.box[d]); + marginDistance2 += pbcMargin*gmx::square(0.5)*norm2(pbc.box[d]); } } } @@ -894,7 +895,8 @@ static bool pullGroupObeysPbcRestrictions(const pull_group_work_t &group, int pullCheckPbcWithinGroups(const pull_t &pull, const rvec *x, - const t_pbc &pbc) + const t_pbc &pbc, + real pbcMargin) { if (pbc.ePBC == epbcNONE) { @@ -924,7 +926,7 @@ int pullCheckPbcWithinGroups(const pull_t &pull, { const pull_group_work_t &group = pull.group[g]; if (group.epgrppbc == epgrppbcREFAT && - !pullGroupObeysPbcRestrictions(group, dimUsed[g], x, pbc, pull.comm.pbcAtomBuffer[g])) + !pullGroupObeysPbcRestrictions(group, dimUsed[g], x, pbc, pull.comm.pbcAtomBuffer[g], pbcMargin)) { return g; }