Convert gmx_mtop_t to C++
[gromacs.git] / src / gromacs / mdlib / constraintrange.cpp
blobabe685a04d34d6f52df220864c5b1c574d4bc192
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 static void constr_recur(const t_blocka *at2con,
56 const t_ilist *ilist, const t_iparams *iparams,
57 gmx_bool bTopB,
58 int at, int depth, int nc, int *path,
59 real r0, real r1, real *r2max,
60 int *count)
62 int ncon1;
63 t_iatom *ia1, *ia2;
64 int c, con, a1;
65 gmx_bool bUse;
66 t_iatom *ia;
67 real len, rn0, rn1;
69 (*count)++;
71 ncon1 = ilist[F_CONSTR].nr/3;
72 ia1 = ilist[F_CONSTR].iatoms;
73 ia2 = ilist[F_CONSTRNC].iatoms;
75 /* Loop over all constraints connected to this atom */
76 for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
78 con = at2con->a[c];
79 /* Do not walk over already used constraints */
80 bUse = TRUE;
81 for (a1 = 0; a1 < depth; a1++)
83 if (con == path[a1])
85 bUse = FALSE;
88 if (bUse)
90 ia = constr_iatomptr(ncon1, ia1, ia2, con);
91 /* Flexible constraints currently have length 0, which is incorrect */
92 if (!bTopB)
94 len = iparams[ia[0]].constr.dA;
96 else
98 len = iparams[ia[0]].constr.dB;
100 /* In the worst case the bond directions alternate */
101 if (nc % 2 == 0)
103 rn0 = r0 + len;
104 rn1 = r1;
106 else
108 rn0 = r0;
109 rn1 = r1 + len;
111 /* Assume angles of 120 degrees between all bonds */
112 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
114 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
115 if (debug)
117 fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
118 for (a1 = 0; a1 < depth; a1++)
120 fprintf(debug, " %d %5.3f",
121 path[a1],
122 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
124 fprintf(debug, " %d %5.3f\n", con, len);
127 /* Limit the number of recursions to 1000*nc,
128 * so a call does not take more than a second,
129 * even for highly connected systems.
131 if (depth + 1 < nc && *count < 1000*nc)
133 if (ia[1] == at)
135 a1 = ia[2];
137 else
139 a1 = ia[1];
141 /* Recursion */
142 path[depth] = con;
143 constr_recur(at2con, ilist, iparams,
144 bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
145 path[depth] = -1;
151 static real constr_r_max_moltype(const gmx_moltype_t *molt,
152 const t_iparams *iparams,
153 const t_inputrec *ir)
155 int natoms, nflexcon, *path, at, count;
157 t_blocka at2con;
158 real r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
160 if (molt->ilist[F_CONSTR].nr == 0 &&
161 molt->ilist[F_CONSTRNC].nr == 0)
163 return 0;
166 natoms = molt->atoms.nr;
168 at2con = make_at2con(0, natoms, molt->ilist, iparams,
169 EI_DYNAMICS(ir->eI), &nflexcon);
170 snew(path, 1+ir->nProjOrder);
171 for (at = 0; at < 1+ir->nProjOrder; at++)
173 path[at] = -1;
176 r2maxA = 0;
177 for (at = 0; at < natoms; at++)
179 r0 = 0;
180 r1 = 0;
182 count = 0;
183 constr_recur(&at2con, molt->ilist, iparams,
184 FALSE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxA, &count);
186 if (ir->efep == efepNO)
188 rmax = sqrt(r2maxA);
190 else
192 r2maxB = 0;
193 for (at = 0; at < natoms; at++)
195 r0 = 0;
196 r1 = 0;
197 count = 0;
198 constr_recur(&at2con, molt->ilist, iparams,
199 TRUE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxB, &count);
201 lam0 = ir->fepvals->init_lambda;
202 if (EI_DYNAMICS(ir->eI))
204 lam0 += ir->init_step*ir->fepvals->delta_lambda;
206 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
207 if (EI_DYNAMICS(ir->eI))
209 lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
210 rmax = std::max(rmax, (1 - lam1)*std::sqrt(r2maxA) + lam1*std::sqrt(r2maxB));
214 done_blocka(&at2con);
215 sfree(path);
217 return rmax;
220 real constr_r_max(FILE *fplog, const gmx_mtop_t *mtop, const t_inputrec *ir)
222 real rmax = 0;
223 for (const gmx_moltype_t &molt : mtop->moltype)
225 rmax = std::max(rmax,
226 constr_r_max_moltype(&molt,
227 mtop->ffparams.iparams, ir));
230 if (fplog)
232 fprintf(fplog, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir->nProjOrder, rmax);
235 return rmax;