Clean up make_at2con
[gromacs.git] / src / gromacs / mdlib / constraintrange.cpp
blob95cf473a8ce2acad822eb309f33424c37ad641dc
1 /*
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! */
38 #include "gmxpre.h"
40 #include "constraintrange.h"
42 #include <cmath>
44 #include <algorithm>
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/real.h"
53 #include "gromacs/utility/smalloc.h"
55 namespace gmx
58 //! Recursing function to help find all adjacent constraints.
59 static void constr_recur(const t_blocka *at2con,
60 const t_ilist *ilist, const t_iparams *iparams,
61 gmx_bool bTopB,
62 int at, int depth, int nc, int *path,
63 real r0, real r1, real *r2max,
64 int *count)
66 int ncon1;
67 t_iatom *ia1, *ia2;
68 int c, con, a1;
69 gmx_bool bUse;
70 t_iatom *ia;
71 real len, rn0, rn1;
73 (*count)++;
75 ncon1 = ilist[F_CONSTR].nr/3;
76 ia1 = ilist[F_CONSTR].iatoms;
77 ia2 = ilist[F_CONSTRNC].iatoms;
79 /* Loop over all constraints connected to this atom */
80 for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
82 con = at2con->a[c];
83 /* Do not walk over already used constraints */
84 bUse = TRUE;
85 for (a1 = 0; a1 < depth; a1++)
87 if (con == path[a1])
89 bUse = FALSE;
92 if (bUse)
94 ia = constr_iatomptr(ncon1, ia1, ia2, con);
95 /* Flexible constraints currently have length 0, which is incorrect */
96 if (!bTopB)
98 len = iparams[ia[0]].constr.dA;
100 else
102 len = iparams[ia[0]].constr.dB;
104 /* In the worst case the bond directions alternate */
105 if (nc % 2 == 0)
107 rn0 = r0 + len;
108 rn1 = r1;
110 else
112 rn0 = r0;
113 rn1 = r1 + len;
115 /* Assume angles of 120 degrees between all bonds */
116 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
118 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
119 if (debug)
121 fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
122 for (a1 = 0; a1 < depth; a1++)
124 fprintf(debug, " %d %5.3f",
125 path[a1],
126 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
128 fprintf(debug, " %d %5.3f\n", con, len);
131 /* Limit the number of recursions to 1000*nc,
132 * so a call does not take more than a second,
133 * even for highly connected systems.
135 if (depth + 1 < nc && *count < 1000*nc)
137 if (ia[1] == at)
139 a1 = ia[2];
141 else
143 a1 = ia[1];
145 /* Recursion */
146 path[depth] = con;
147 constr_recur(at2con, ilist, iparams,
148 bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
149 path[depth] = -1;
155 //! Find the interaction radius needed for constraints for this molecule type.
156 static real constr_r_max_moltype(const gmx_moltype_t *molt,
157 const t_iparams *iparams,
158 const t_inputrec *ir)
160 int natoms, *path, at, count;
162 t_blocka at2con;
163 real r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
165 if (molt->ilist[F_CONSTR].nr == 0 &&
166 molt->ilist[F_CONSTRNC].nr == 0)
168 return 0;
171 natoms = molt->atoms.nr;
173 at2con = make_at2con(natoms, molt->ilist, iparams,
174 flexibleConstraintTreatment(EI_DYNAMICS(ir->eI)));
175 snew(path, 1+ir->nProjOrder);
176 for (at = 0; at < 1+ir->nProjOrder; at++)
178 path[at] = -1;
181 r2maxA = 0;
182 for (at = 0; at < natoms; at++)
184 r0 = 0;
185 r1 = 0;
187 count = 0;
188 constr_recur(&at2con, molt->ilist, iparams,
189 FALSE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxA, &count);
191 if (ir->efep == efepNO)
193 rmax = sqrt(r2maxA);
195 else
197 r2maxB = 0;
198 for (at = 0; at < natoms; at++)
200 r0 = 0;
201 r1 = 0;
202 count = 0;
203 constr_recur(&at2con, molt->ilist, iparams,
204 TRUE, at, 0, 1+ir->nProjOrder, path, r0, 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));
219 done_blocka(&at2con);
220 sfree(path);
222 return rmax;
225 real constr_r_max(FILE *fplog, const gmx_mtop_t *mtop, const t_inputrec *ir)
227 real rmax = 0;
228 for (const gmx_moltype_t &molt : mtop->moltype)
230 rmax = std::max(rmax,
231 constr_r_max_moltype(&molt,
232 mtop->ffparams.iparams, ir));
235 if (fplog)
237 fprintf(fplog, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir->nProjOrder, rmax);
240 return rmax;
243 } // namespace