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