Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / topshake.c
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, 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 "topshake.h"
41
42 #include <ctype.h>
43 #include <math.h>
44
45 #include "gromacs/gmxpreprocess/readir.h"
46 #include "gromacs/gmxpreprocess/topdirs.h"
47 #include "gromacs/gmxpreprocess/toppush.h"
48 #include "gromacs/gmxpreprocess/toputil.h"
49 #include "gromacs/legacyheaders/typedefs.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
53
54 static void copy_bond (t_params *pr, int to, int from)
55 /* copies an entry in a bond list to another position.
56  * does no allocing or freeing of memory
57  */
58 {
59     /*memcpy((char*) &(pr->param[to]),(char*) &(pr->param[from]),
60        (size_t)sizeof(pr->param[from]));*/
61     int i;
62
63     if (to != from)
64     {
65         range_check(to, 0, pr->nr);
66         range_check(from, 0, pr->nr);
67         for (i = 0; (i < MAXATOMLIST); i++)
68         {
69             pr->param[to].a[i] = pr->param[from].a[i];
70         }
71         for (i = 0; (i < MAXFORCEPARAM); i++)
72         {
73             pr->param[to].c[i] = pr->param[from].c[i];
74         }
75         for (i = 0; (i < MAXSLEN); i++)
76         {
77             pr->param[to].s[i] = pr->param[from].s[i];
78         }
79     }
80 }
81
82 static int count_hydrogens (char ***atomname, int nra, atom_id a[])
83 {
84     int  i, nh;
85
86     if (!atomname)
87     {
88         gmx_fatal(FARGS, "Cannot call count_hydrogens with no atomname (%s %d)",
89                   __FILE__, __LINE__);
90     }
91
92     nh = 0;
93     for (i = 0; (i < nra); i++)
94     {
95         if (toupper(**(atomname[a[i]])) == 'H')
96         {
97             nh++;
98         }
99     }
100
101     return nh;
102 }
103
104 void make_shake (t_params plist[], t_atoms *atoms, int nshake)
105 {
106     char          ***info = atoms->atomname;
107     t_params        *pr;
108     t_params        *bonds;
109     t_param          p, *bond, *ang;
110     real             b_ij, b_jk;
111     int              nb, b, i, j, ftype, ftype_a;
112     gmx_bool         bFound;
113
114     if (nshake != eshNONE)
115     {
116         switch (nshake)
117         {
118             case eshHBONDS:
119                 printf("turning H bonds into constraints...\n");
120                 break;
121             case eshALLBONDS:
122                 printf("turning all bonds into constraints...\n");
123                 break;
124             case eshHANGLES:
125                 printf("turning all bonds and H angles into constraints...\n");
126                 break;
127             case eshALLANGLES:
128                 printf("turning all bonds and angles into constraints...\n");
129                 break;
130             default:
131                 gmx_fatal(FARGS, "Invalid option for make_shake (%d)", nshake);
132         }
133
134         if ((nshake == eshHANGLES) || (nshake == eshALLANGLES))
135         {
136             /* Add all the angles with hydrogens to the shake list
137              * and remove them from the bond list
138              */
139             for (ftype = 0; (ftype < F_NRE); ftype++)
140             {
141                 if (interaction_function[ftype].flags & IF_BTYPE)
142                 {
143                     bonds = &(plist[ftype]);
144
145                     for (ftype_a = 0; (bonds->nr > 0 && ftype_a < F_NRE); ftype_a++)
146                     {
147                         if (interaction_function[ftype_a].flags & IF_ATYPE)
148                         {
149                             pr = &(plist[ftype_a]);
150
151                             for (i = 0; (i < pr->nr); )
152                             {
153                                 int numhydrogens;
154
155                                 ang = &(pr->param[i]);
156 #ifdef DEBUG
157                                 printf("Angle: %d-%d-%d\n", ang->AI, ang->AJ, ang->AK);
158 #endif
159                                 numhydrogens = count_hydrogens(info, 3, ang->a);
160                                 if ((nshake == eshALLANGLES) ||
161                                     (numhydrogens > 1) ||
162                                     (numhydrogens == 1 && toupper(**(info[ang->a[1]])) == 'O'))
163                                 {
164                                     /* Can only add hydrogen angle shake, if the two bonds
165                                      * are constrained.
166                                      * append this angle to the shake list
167                                      */
168                                     p.AI = ang->AI;
169                                     p.AJ = ang->AK;
170
171                                     /* Calculate length of constraint */
172                                     bFound = FALSE;
173                                     b_ij   = b_jk = 0.0;
174                                     for (j = 0; !bFound && (j < bonds->nr); j++)
175                                     {
176                                         bond = &(bonds->param[j]);
177                                         if (((bond->AI == ang->AI) &&
178                                              (bond->AJ == ang->AJ)) ||
179                                             ((bond->AI == ang->AJ) &&
180                                              (bond->AJ == ang->AI)))
181                                         {
182                                             b_ij = bond->C0;
183                                         }
184                                         if (((bond->AI == ang->AK) &&
185                                              (bond->AJ == ang->AJ)) ||
186                                             ((bond->AI == ang->AJ) &&
187                                              (bond->AJ == ang->AK)))
188                                         {
189                                             b_jk = bond->C0;
190                                         }
191                                         bFound = (b_ij != 0.0) && (b_jk != 0.0);
192                                     }
193                                     if (bFound)
194                                     {
195                                         /* apply law of cosines */
196                                         p.C0 = sqrt( b_ij*b_ij + b_jk*b_jk -
197                                                      2.0*b_ij*b_jk*cos(DEG2RAD*ang->C0) );
198                                         p.C1 = p.C0;
199 #ifdef DEBUG
200                                         printf("p: %d, q: %d, dist: %12.5e\n", p.AI, p.AJ, p.C0);
201 #endif
202                                         add_param_to_list (&(plist[F_CONSTR]), &p);
203                                         /* move the last bond to this position */
204                                         copy_bond (pr, i, pr->nr-1);
205                                         /* should free memory here!! */
206                                         pr->nr--;
207                                     }
208                                 }
209                                 else
210                                 {
211                                     i++;
212                                 }
213                             }
214                         } /* if IF_ATYPE */
215                     }     /* for ftype_A */
216                 }         /* if IF_BTYPE */
217             }             /* for ftype */
218         }                 /* if shake angles */
219
220         /* Add all the bonds with hydrogens to the shake list
221          * and remove them from the bond list
222          */
223         for (ftype = 0; (ftype < F_NRE); ftype++)
224         {
225             if (interaction_function[ftype].flags & IF_BTYPE)
226             {
227                 pr = &(plist[ftype]);
228                 j  = 0;
229                 for (i = 0; i < pr->nr; i++)
230                 {
231                     if ( (nshake != eshHBONDS) ||
232                          (count_hydrogens (info, 2, pr->param[i].a) > 0) )
233                     {
234                         /* append this bond to the shake list */
235                         p.AI = pr->param[i].AI;
236                         p.AJ = pr->param[i].AJ;
237                         p.C0 = pr->param[i].C0;
238                         p.C1 = pr->param[i].C2;
239                         add_param_to_list (&(plist[F_CONSTR]), &p);
240                     }
241                     else
242                     {
243                         copy_bond(pr, j++, i);
244                     }
245                 }
246                 pr->nr = j;
247             }
248         }
249     }
250 }