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