Use ListOfLists for atom to constraint mapping
[alexxy/gromacs.git] / src / gromacs / mdlib / constraintrange.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2019, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "constraintrange.h"
41
42 #include <cmath>
43
44 #include <algorithm>
45
46 #include "gromacs/mdlib/constr.h"
47 #include "gromacs/mdtypes/inputrec.h"
48 #include "gromacs/topology/mtop_util.h"
49 #include "gromacs/utility/basedefinitions.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/listoflists.h"
52 #include "gromacs/utility/logger.h"
53 #include "gromacs/utility/real.h"
54
55 namespace gmx
56 {
57
58 //! Recursing function to help find all adjacent constraints.
59 static void constr_recur(const ListOfLists<int>&        at2con,
60                          const InteractionLists&        ilist,
61                          gmx::ArrayRef<const t_iparams> iparams,
62                          gmx_bool                       bTopB,
63                          int                            at,
64                          int                            depth,
65                          int                            nc,
66                          ArrayRef<int>                  path,
67                          real                           r0,
68                          real                           r1,
69                          real*                          r2max,
70                          int*                           count)
71 {
72     gmx_bool bUse;
73     real     len, rn0, rn1;
74
75     (*count)++;
76
77     gmx::ArrayRef<const int> ia1 = ilist[F_CONSTR].iatoms;
78     gmx::ArrayRef<const int> ia2 = ilist[F_CONSTRNC].iatoms;
79
80     /* Loop over all constraints connected to this atom */
81     for (const int con : at2con[at])
82     {
83         /* Do not walk over already used constraints */
84         bUse = TRUE;
85         for (int a1 = 0; a1 < depth; a1++)
86         {
87             if (con == path[a1])
88             {
89                 bUse = FALSE;
90             }
91         }
92         if (bUse)
93         {
94             const int* ia = constr_iatomptr(ia1, ia2, con);
95             /* Flexible constraints currently have length 0, which is incorrect */
96             if (!bTopB)
97             {
98                 len = iparams[ia[0]].constr.dA;
99             }
100             else
101             {
102                 len = iparams[ia[0]].constr.dB;
103             }
104             /* In the worst case the bond directions alternate */
105             if (nc % 2 == 0)
106             {
107                 rn0 = r0 + len;
108                 rn1 = r1;
109             }
110             else
111             {
112                 rn0 = r0;
113                 rn1 = r1 + len;
114             }
115             /* Assume angles of 120 degrees between all bonds */
116             if (rn0 * rn0 + rn1 * rn1 + rn0 * rn1 > *r2max)
117             {
118                 *r2max = rn0 * rn0 + rn1 * rn1 + r0 * rn1;
119                 if (debug)
120                 {
121                     fprintf(debug,
122                             "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0,
123                             rn1, sqrt(*r2max));
124                     for (int a1 = 0; a1 < depth; a1++)
125                     {
126                         fprintf(debug, " %d %5.3f", path[a1],
127                                 iparams[constr_iatomptr(ia1, ia2, con)[0]].constr.dA);
128                     }
129                     fprintf(debug, " %d %5.3f\n", con, len);
130                 }
131             }
132             /* Limit the number of recursions to 1000*nc,
133              * so a call does not take more than a second,
134              * even for highly connected systems.
135              */
136             if (depth + 1 < nc && *count < 1000 * nc)
137             {
138                 int a1;
139                 if (ia[1] == at)
140                 {
141                     a1 = ia[2];
142                 }
143                 else
144                 {
145                     a1 = ia[1];
146                 }
147                 /* Recursion */
148                 path[depth] = con;
149                 constr_recur(at2con, ilist, iparams, bTopB, a1, depth + 1, nc, path, rn0, rn1, r2max, count);
150                 path[depth] = -1;
151             }
152         }
153     }
154 }
155
156 //! Find the interaction radius needed for constraints for this molecule type.
157 static real constr_r_max_moltype(const gmx_moltype_t*           molt,
158                                  gmx::ArrayRef<const t_iparams> iparams,
159                                  const t_inputrec*              ir)
160 {
161     int natoms, at, count;
162
163     real r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
164
165     if (molt->ilist[F_CONSTR].size() == 0 && molt->ilist[F_CONSTRNC].size() == 0)
166     {
167         return 0;
168     }
169
170     natoms = molt->atoms.nr;
171
172     const ListOfLists<int> at2con =
173             make_at2con(*molt, iparams, flexibleConstraintTreatment(EI_DYNAMICS(ir->eI)));
174     std::vector<int> path(1 + ir->nProjOrder);
175     for (at = 0; at < 1 + ir->nProjOrder; at++)
176     {
177         path[at] = -1;
178     }
179
180     r2maxA = 0;
181     for (at = 0; at < natoms; at++)
182     {
183         r0 = 0;
184         r1 = 0;
185
186         count = 0;
187         constr_recur(at2con, molt->ilist, iparams, FALSE, at, 0, 1 + ir->nProjOrder, path, r0, r1,
188                      &r2maxA, &count);
189     }
190     if (ir->efep == efepNO)
191     {
192         rmax = sqrt(r2maxA);
193     }
194     else
195     {
196         r2maxB = 0;
197         for (at = 0; at < natoms; at++)
198         {
199             r0    = 0;
200             r1    = 0;
201             count = 0;
202             constr_recur(at2con, molt->ilist, iparams, TRUE, at, 0, 1 + ir->nProjOrder, path, r0,
203                          r1, &r2maxB, &count);
204         }
205         lam0 = ir->fepvals->init_lambda;
206         if (EI_DYNAMICS(ir->eI))
207         {
208             lam0 += ir->init_step * ir->fepvals->delta_lambda;
209         }
210         rmax = (1 - lam0) * sqrt(r2maxA) + lam0 * sqrt(r2maxB);
211         if (EI_DYNAMICS(ir->eI))
212         {
213             lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps) * ir->fepvals->delta_lambda;
214             rmax = std::max(rmax, (1 - lam1) * std::sqrt(r2maxA) + lam1 * std::sqrt(r2maxB));
215         }
216     }
217
218     return rmax;
219 }
220
221 real constr_r_max(const MDLogger& mdlog, const gmx_mtop_t* mtop, const t_inputrec* ir)
222 {
223     real rmax = 0;
224     for (const gmx_moltype_t& molt : mtop->moltype)
225     {
226         rmax = std::max(rmax, constr_r_max_moltype(&molt, mtop->ffparams.iparams, ir));
227     }
228
229     GMX_LOG(mdlog.info)
230             .appendTextFormatted(
231                     "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm",
232                     1 + ir->nProjOrder, rmax);
233
234     return rmax;
235 }
236
237 } // namespace gmx