Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxana / nrama.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 #include "gmxpre.h"
38
39 #include "nrama.h"
40
41 #include <math.h>
42 #include <stdlib.h>
43
44 #include "gromacs/bonded/bonded.h"
45 #include "gromacs/pbcutil/rmpbc.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/futil.h"
48 #include "gromacs/utility/smalloc.h"
49
50 static const char *pp_pat[] = { "C", "N", "CA", "C", "N" };
51 #define NPP (sizeof(pp_pat)/sizeof(pp_pat[0]))
52
53 static int d_comp(const void *a, const void *b)
54 {
55     t_dih *da, *db;
56
57     da = (t_dih *)a;
58     db = (t_dih *)b;
59
60     if (da->ai[1] < db->ai[1])
61     {
62         return -1;
63     }
64     else if (da->ai[1] == db->ai[1])
65     {
66         return (da->ai[2] - db->ai[2]);
67     }
68     else
69     {
70         return 1;
71     }
72 }
73
74
75 static void calc_dihs(t_xrama *xr)
76 {
77     int          i, t1, t2, t3;
78     rvec         r_ij, r_kj, r_kl, m, n;
79     real         sign;
80     t_dih       *dd;
81     gmx_rmpbc_t  gpbc = NULL;
82
83     gpbc = gmx_rmpbc_init(xr->idef, xr->ePBC, xr->natoms);
84     gmx_rmpbc(gpbc, xr->natoms, xr->box, xr->x);
85     gmx_rmpbc_done(gpbc);
86
87     for (i = 0; (i < xr->ndih); i++)
88     {
89         dd      = &(xr->dih[i]);
90         dd->ang = dih_angle(xr->x[dd->ai[0]], xr->x[dd->ai[1]],
91                             xr->x[dd->ai[2]], xr->x[dd->ai[3]],
92                             NULL,
93                             r_ij, r_kj, r_kl, m, n, &sign, &t1, &t2, &t3);
94     }
95 }
96
97 gmx_bool new_data(t_xrama *xr)
98 {
99     if (!read_next_x(xr->oenv, xr->traj, &xr->t, xr->x, xr->box))
100     {
101         return FALSE;
102     }
103
104     calc_dihs(xr);
105
106     return TRUE;
107 }
108
109 static int find_atom(const char *find, char ***names, int start, int nr)
110 {
111     int i;
112
113     for (i = start; (i < nr); i++)
114     {
115         if (strcmp(find, *names[i]) == 0)
116         {
117             return i;
118         }
119     }
120     return -1;
121 }
122
123 static void add_xr(t_xrama *xr, int ff[5], t_atoms *atoms)
124 {
125     char buf[12];
126     int  i;
127
128     srenew(xr->dih, xr->ndih+2);
129     for (i = 0; (i < 4); i++)
130     {
131         xr->dih[xr->ndih].ai[i] = ff[i];
132     }
133     for (i = 0; (i < 4); i++)
134     {
135         xr->dih[xr->ndih+1].ai[i] = ff[i+1];
136     }
137     xr->ndih += 2;
138
139     srenew(xr->pp, xr->npp+1);
140     xr->pp[xr->npp].iphi  = xr->ndih-2;
141     xr->pp[xr->npp].ipsi  = xr->ndih-1;
142     xr->pp[xr->npp].bShow = FALSE;
143     sprintf(buf, "%s-%d", *atoms->resinfo[atoms->atom[ff[1]].resind].name,
144             atoms->resinfo[atoms->atom[ff[1]].resind].nr);
145     xr->pp[xr->npp].label = gmx_strdup(buf);
146     xr->npp++;
147 }
148
149 static void get_dih(t_xrama *xr, t_atoms *atoms)
150 {
151     int    found, ff[NPP];
152     int    i;
153     size_t j;
154
155     for (i = 0; (i < atoms->nr); )
156     {
157         found = i;
158         for (j = 0; (j < NPP); j++)
159         {
160             if ((ff[j] = find_atom(pp_pat[j], atoms->atomname, found, atoms->nr)) == -1)
161             {
162                 break;
163             }
164             found = ff[j]+1;
165         }
166         if (j != NPP)
167         {
168             break;
169         }
170         add_xr(xr, ff, atoms);
171         i = ff[0]+1;
172     }
173     fprintf(stderr, "Found %d phi-psi combinations\n", xr->npp);
174 }
175
176 static int search_ff(int thisff[NPP], int ndih, int **ff)
177 {
178     int      j, k;
179     gmx_bool bFound = FALSE;
180
181     for (j = 0; (j < ndih); j++)
182     {
183         bFound = TRUE;
184         for (k = 1; (k <= 3); k++)
185         {
186             bFound = bFound && (thisff[k] == ff[j][k]);
187         }
188         if (bFound)
189         {
190             if (thisff[0] == -1)
191             {
192                 ff[j][4] = thisff[4];
193             }
194             else
195             {
196                 ff[j][0] = thisff[0];
197             }
198             break;
199         }
200     }
201     if (!bFound)
202     {
203         for (k = 0; (k < 5); k++)
204         {
205             ff[ndih][k] = thisff[k];
206         }
207         ndih++;
208     }
209     return ndih;
210 }
211
212 static void min_max(t_xrama *xr)
213 {
214     int ai, i, j;
215
216     xr->amin = xr->natoms;
217     xr->amax = 0;
218     for (i = 0; (i < xr->ndih); i++)
219     {
220         for (j = 0; (j < 4); j++)
221         {
222             ai = xr->dih[i].ai[j];
223             if (ai < xr->amin)
224             {
225                 xr->amin = ai;
226             }
227             else if (ai > xr->amax)
228             {
229                 xr->amax = ai;
230             }
231         }
232     }
233 }
234
235 static void get_dih_props(t_xrama *xr, t_idef *idef, int mult)
236 {
237     int      i, ft, ftype, nra;
238     t_iatom *ia;
239     t_dih   *dd, key;
240
241     ia = idef->il[F_PDIHS].iatoms;
242     for (i = 0; (i < idef->il[F_PDIHS].nr); )
243     {
244         ft    = ia[0];
245         ftype = idef->functype[ft];
246         nra   = interaction_function[ftype].nratoms;
247
248         if (ftype != F_PDIHS)
249         {
250             gmx_incons("ftype is not a dihedral");
251         }
252
253         key.ai[1] = ia[2];
254         key.ai[2] = ia[3];
255         if ((dd = (t_dih *)bsearch(&key, xr->dih, xr->ndih, (size_t)sizeof(key), d_comp))
256             != NULL)
257         {
258             dd->mult = idef->iparams[ft].pdihs.mult;
259             dd->phi0 = idef->iparams[ft].pdihs.phiA;
260         }
261
262         i  += nra+1;
263         ia += nra+1;
264     }
265     /* Fill in defaults for values not in the topology */
266     for (i = 0; (i < xr->ndih); i++)
267     {
268         if (xr->dih[i].mult == 0)
269         {
270             fprintf(stderr,
271                     "Dihedral around %d,%d not found in topology. Using mult=%d\n",
272                     xr->dih[i].ai[1], xr->dih[i].ai[2], mult);
273             xr->dih[i].mult = mult;
274             xr->dih[i].phi0 = 180;
275         }
276     }
277 }
278
279
280
281 t_topology *init_rama(const output_env_t oenv, const char *infile,
282                       const char *topfile, t_xrama *xr, int mult)
283 {
284     t_topology *top;
285     int         ePBC;
286     real        t;
287
288     top = read_top(topfile, &xr->ePBC);
289
290     /*get_dih2(xr,top->idef.functype,&(top->idef.bondeds),&(top->atoms));*/
291     get_dih(xr, &(top->atoms));
292     get_dih_props(xr, &(top->idef), mult);
293     xr->natoms = read_first_x(oenv, &xr->traj, infile, &t, &(xr->x), xr->box);
294     xr->idef   = &(top->idef);
295     xr->oenv   = oenv;
296
297     min_max(xr);
298     calc_dihs(xr);
299
300     return top;
301 }