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