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) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source 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.
38 /* This file is completely threadsafe - keep it that way! */
41 #include "constraintrange.h"
47 #include "gromacs/mdlib/constr.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/topology/mtop_util.h"
50 #include "gromacs/utility/basedefinitions.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/listoflists.h"
53 #include "gromacs/utility/logger.h"
54 #include "gromacs/utility/real.h"
59 //! Recursing function to help find all adjacent constraints.
60 static void constr_recur(const ListOfLists
<int>& at2con
,
61 const InteractionLists
& ilist
,
62 gmx::ArrayRef
<const t_iparams
> iparams
,
78 gmx::ArrayRef
<const int> ia1
= ilist
[F_CONSTR
].iatoms
;
79 gmx::ArrayRef
<const int> ia2
= ilist
[F_CONSTRNC
].iatoms
;
81 /* Loop over all constraints connected to this atom */
82 for (const int con
: at2con
[at
])
84 /* Do not walk over already used constraints */
86 for (int a1
= 0; a1
< depth
; a1
++)
95 const int* ia
= constr_iatomptr(ia1
, ia2
, con
);
96 /* Flexible constraints currently have length 0, which is incorrect */
99 len
= iparams
[ia
[0]].constr
.dA
;
103 len
= iparams
[ia
[0]].constr
.dB
;
105 /* In the worst case the bond directions alternate */
116 /* Assume angles of 120 degrees between all bonds */
117 if (rn0
* rn0
+ rn1
* rn1
+ rn0
* rn1
> *r2max
)
119 *r2max
= rn0
* rn0
+ rn1
* rn1
+ r0
* rn1
;
123 "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0
,
125 for (int a1
= 0; a1
< depth
; a1
++)
127 fprintf(debug
, " %d %5.3f", path
[a1
],
128 iparams
[constr_iatomptr(ia1
, ia2
, con
)[0]].constr
.dA
);
130 fprintf(debug
, " %d %5.3f\n", con
, len
);
133 /* Limit the number of recursions to 1000*nc,
134 * so a call does not take more than a second,
135 * even for highly connected systems.
137 if (depth
+ 1 < nc
&& *count
< 1000 * nc
)
150 constr_recur(at2con
, ilist
, iparams
, bTopB
, a1
, depth
+ 1, nc
, path
, rn0
, rn1
, r2max
, count
);
157 //! Find the interaction radius needed for constraints for this molecule type.
158 static real
constr_r_max_moltype(const gmx_moltype_t
* molt
,
159 gmx::ArrayRef
<const t_iparams
> iparams
,
160 const t_inputrec
* ir
)
162 int natoms
, at
, count
;
164 real r0
, r1
, r2maxA
, r2maxB
, rmax
, lam0
, lam1
;
166 if (molt
->ilist
[F_CONSTR
].empty() && molt
->ilist
[F_CONSTRNC
].empty())
171 natoms
= molt
->atoms
.nr
;
173 const ListOfLists
<int> at2con
=
174 make_at2con(*molt
, iparams
, flexibleConstraintTreatment(EI_DYNAMICS(ir
->eI
)));
175 std::vector
<int> path(1 + ir
->nProjOrder
);
176 for (at
= 0; at
< 1 + ir
->nProjOrder
; at
++)
182 for (at
= 0; at
< natoms
; at
++)
188 constr_recur(at2con
, molt
->ilist
, iparams
, FALSE
, at
, 0, 1 + ir
->nProjOrder
, path
, r0
, r1
,
191 if (ir
->efep
== efepNO
)
198 for (at
= 0; at
< natoms
; at
++)
203 constr_recur(at2con
, molt
->ilist
, iparams
, TRUE
, at
, 0, 1 + ir
->nProjOrder
, path
, r0
,
204 r1
, &r2maxB
, &count
);
206 lam0
= ir
->fepvals
->init_lambda
;
207 if (EI_DYNAMICS(ir
->eI
))
209 lam0
+= ir
->init_step
* ir
->fepvals
->delta_lambda
;
211 rmax
= (1 - lam0
) * sqrt(r2maxA
) + lam0
* sqrt(r2maxB
);
212 if (EI_DYNAMICS(ir
->eI
))
214 lam1
= ir
->fepvals
->init_lambda
+ (ir
->init_step
+ ir
->nsteps
) * ir
->fepvals
->delta_lambda
;
215 rmax
= std::max(rmax
, (1 - lam1
) * std::sqrt(r2maxA
) + lam1
* std::sqrt(r2maxB
));
222 real
constr_r_max(const MDLogger
& mdlog
, const gmx_mtop_t
* mtop
, const t_inputrec
* ir
)
225 for (const gmx_moltype_t
& molt
: mtop
->moltype
)
227 rmax
= std::max(rmax
, constr_r_max_moltype(&molt
, mtop
->ffparams
.iparams
, ir
));
231 .appendTextFormatted(
232 "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm",
233 1 + ir
->nProjOrder
, rmax
);