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,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /* This file is completely threadsafe - keep it that way! */
40 #include "constraintrange.h"
46 #include "gromacs/mdlib/constr.h"
47 #include "gromacs/mdtypes/inputrec.h"
48 #include "gromacs/topology/block.h"
49 #include "gromacs/topology/mtop_util.h"
50 #include "gromacs/utility/basedefinitions.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/logger.h"
53 #include "gromacs/utility/real.h"
54 #include "gromacs/utility/smalloc.h"
59 //! Recursing function to help find all adjacent constraints.
60 static void constr_recur(const t_blocka
*at2con
,
61 const InteractionLists
&ilist
,
62 gmx::ArrayRef
<const t_iparams
> iparams
,
64 int at
, int depth
, int nc
, int *path
,
65 real r0
, real r1
, real
*r2max
,
74 gmx::ArrayRef
<const int> ia1
= ilist
[F_CONSTR
].iatoms
;
75 gmx::ArrayRef
<const int> ia2
= ilist
[F_CONSTRNC
].iatoms
;
77 /* Loop over all constraints connected to this atom */
78 for (c
= at2con
->index
[at
]; c
< at2con
->index
[at
+1]; c
++)
81 /* Do not walk over already used constraints */
83 for (a1
= 0; a1
< depth
; a1
++)
92 const int *ia
= constr_iatomptr(ia1
, ia2
, con
);
93 /* Flexible constraints currently have length 0, which is incorrect */
96 len
= iparams
[ia
[0]].constr
.dA
;
100 len
= iparams
[ia
[0]].constr
.dB
;
102 /* In the worst case the bond directions alternate */
113 /* Assume angles of 120 degrees between all bonds */
114 if (rn0
*rn0
+ rn1
*rn1
+ rn0
*rn1
> *r2max
)
116 *r2max
= rn0
*rn0
+ rn1
*rn1
+ r0
*rn1
;
119 fprintf(debug
, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0
, rn1
, sqrt(*r2max
));
120 for (a1
= 0; a1
< depth
; a1
++)
122 fprintf(debug
, " %d %5.3f",
124 iparams
[constr_iatomptr(ia1
, ia2
, con
)[0]].constr
.dA
);
126 fprintf(debug
, " %d %5.3f\n", con
, len
);
129 /* Limit the number of recursions to 1000*nc,
130 * so a call does not take more than a second,
131 * even for highly connected systems.
133 if (depth
+ 1 < nc
&& *count
< 1000*nc
)
145 constr_recur(at2con
, ilist
, iparams
,
146 bTopB
, a1
, depth
+1, nc
, path
, rn0
, rn1
, r2max
, count
);
153 //! Find the interaction radius needed for constraints for this molecule type.
154 static real
constr_r_max_moltype(const gmx_moltype_t
*molt
,
155 gmx::ArrayRef
<const t_iparams
> iparams
,
156 const t_inputrec
*ir
)
158 int natoms
, *path
, at
, count
;
161 real r0
, r1
, r2maxA
, r2maxB
, rmax
, lam0
, lam1
;
163 if (molt
->ilist
[F_CONSTR
].size() == 0 &&
164 molt
->ilist
[F_CONSTRNC
].size() == 0)
169 natoms
= molt
->atoms
.nr
;
171 at2con
= make_at2con(*molt
, iparams
,
172 flexibleConstraintTreatment(EI_DYNAMICS(ir
->eI
)));
173 snew(path
, 1+ir
->nProjOrder
);
174 for (at
= 0; at
< 1+ir
->nProjOrder
; at
++)
180 for (at
= 0; at
< natoms
; at
++)
186 constr_recur(&at2con
, molt
->ilist
, iparams
,
187 FALSE
, at
, 0, 1+ir
->nProjOrder
, path
, r0
, r1
, &r2maxA
, &count
);
189 if (ir
->efep
== efepNO
)
196 for (at
= 0; at
< natoms
; at
++)
201 constr_recur(&at2con
, molt
->ilist
, iparams
,
202 TRUE
, at
, 0, 1+ir
->nProjOrder
, path
, r0
, r1
, &r2maxB
, &count
);
204 lam0
= ir
->fepvals
->init_lambda
;
205 if (EI_DYNAMICS(ir
->eI
))
207 lam0
+= ir
->init_step
*ir
->fepvals
->delta_lambda
;
209 rmax
= (1 - lam0
)*sqrt(r2maxA
) + lam0
*sqrt(r2maxB
);
210 if (EI_DYNAMICS(ir
->eI
))
212 lam1
= ir
->fepvals
->init_lambda
+ (ir
->init_step
+ ir
->nsteps
)*ir
->fepvals
->delta_lambda
;
213 rmax
= std::max(rmax
, (1 - lam1
)*std::sqrt(r2maxA
) + lam1
*std::sqrt(r2maxB
));
217 done_blocka(&at2con
);
223 real
constr_r_max(const MDLogger
&mdlog
, const gmx_mtop_t
*mtop
, const t_inputrec
*ir
)
226 for (const gmx_moltype_t
&molt
: mtop
->moltype
)
228 rmax
= std::max(rmax
,
229 constr_r_max_moltype(&molt
,
230 mtop
->ffparams
.iparams
, ir
));
233 GMX_LOG(mdlog
.info
).appendTextFormatted(
234 "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm",
235 1+ir
->nProjOrder
, rmax
);