Merge branch 'origin/release-2020' into master
[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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "nrama.h"
41
42 #include <cstdlib>
43 #include <cstring>
44
45 #include <algorithm>
46
47 #include "gromacs/listed_forces/bonded.h"
48 #include "gromacs/pbcutil/rmpbc.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/futil.h"
52 #include "gromacs/utility/smalloc.h"
53
54 static const char* pp_pat[] = { "C", "N", "CA", "C", "N" };
55 #define NPP (sizeof(pp_pat) / sizeof(pp_pat[0]))
56
57 static bool d_comp(const t_dih& a, const t_dih& b)
58 {
59     if (a.ai[1] < b.ai[1])
60     {
61         return true;
62     }
63     else if (a.ai[1] == b.ai[1])
64     {
65         return a.ai[2] < b.ai[2];
66     }
67     else
68     {
69         return false;
70     }
71 }
72
73
74 static void calc_dihs(t_xrama* xr)
75 {
76     int         i, t1, t2, t3;
77     rvec        r_ij, r_kj, r_kl, m, n;
78     t_dih*      dd;
79     gmx_rmpbc_t gpbc = nullptr;
80
81     gpbc = gmx_rmpbc_init(xr->idef, xr->pbcType, 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]], xr->x[dd->ai[2]], xr->x[dd->ai[3]],
89                             nullptr, r_ij, r_kj, r_kl, m, n, &t1, &t2, &t3);
90     }
91 }
92
93 gmx_bool new_data(t_xrama* xr)
94 {
95     if (!read_next_x(xr->oenv, xr->traj, &xr->t, xr->x, xr->box))
96     {
97         return FALSE;
98     }
99
100     calc_dihs(xr);
101
102     return TRUE;
103 }
104
105 static int find_atom(const char* find, char*** names, int start, int nr)
106 {
107     int i;
108
109     for (i = start; (i < nr); i++)
110     {
111         if (std::strcmp(find, *names[i]) == 0)
112         {
113             return i;
114         }
115     }
116     return -1;
117 }
118
119 static void add_xr(t_xrama* xr, const int ff[5], const t_atoms* atoms)
120 {
121     char buf[12];
122     int  i;
123
124     srenew(xr->dih, xr->ndih + 2LL);
125     for (i = 0; (i < 4); i++)
126     {
127         xr->dih[xr->ndih].ai[i] = ff[i];
128     }
129     for (i = 0; (i < 4); i++)
130     {
131         xr->dih[xr->ndih + 1].ai[i] = ff[i + 1];
132     }
133     xr->ndih += 2;
134
135     srenew(xr->pp, xr->npp + 1LL);
136     xr->pp[xr->npp].iphi  = xr->ndih - 2;
137     xr->pp[xr->npp].ipsi  = xr->ndih - 1;
138     xr->pp[xr->npp].bShow = FALSE;
139     sprintf(buf, "%s-%d", *atoms->resinfo[atoms->atom[ff[1]].resind].name,
140             atoms->resinfo[atoms->atom[ff[1]].resind].nr);
141     xr->pp[xr->npp].label = gmx_strdup(buf);
142     xr->npp++;
143 }
144
145 static void get_dih(t_xrama* xr, const t_atoms* atoms)
146 {
147     int    found, ff[NPP];
148     int    i;
149     size_t j;
150
151     for (i = 0; (i < atoms->nr);)
152     {
153         found = i;
154         for (j = 0; (j < NPP); j++)
155         {
156             if ((ff[j] = find_atom(pp_pat[j], atoms->atomname, found, atoms->nr)) == -1)
157             {
158                 break;
159             }
160             found = ff[j] + 1;
161         }
162         if (j != NPP)
163         {
164             break;
165         }
166         add_xr(xr, ff, atoms);
167         i = ff[0] + 1;
168     }
169     fprintf(stderr, "Found %d phi-psi combinations\n", xr->npp);
170 }
171
172 static void min_max(t_xrama* xr)
173 {
174     int ai, i, j;
175
176     xr->amin = xr->natoms;
177     xr->amax = 0;
178     for (i = 0; (i < xr->ndih); i++)
179     {
180         for (j = 0; (j < 4); j++)
181         {
182             MSVC_DIAGNOSTIC_IGNORE(28182) // false positive in 2019 (16.5.4)
183             ai = xr->dih[i].ai[j];
184             MSVC_DIAGNOSTIC_RESET
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 + 1LL;
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, "Dihedral around %d,%d not found in topology. Using mult=%d\n",
233                     xr->dih[i].ai[1], xr->dih[i].ai[2], mult);
234             xr->dih[i].mult = mult;
235             xr->dih[i].phi0 = 180;
236         }
237     }
238 }
239
240
241 t_topology* init_rama(gmx_output_env_t* oenv, const char* infile, const char* topfile, t_xrama* xr, int mult)
242 {
243     t_topology* top;
244     real        t;
245
246     top = read_top(topfile, &xr->pbcType);
247
248     /*get_dih2(xr,top->idef.functype,&(top->idef.bondeds),&(top->atoms));*/
249     get_dih(xr, &(top->atoms));
250     get_dih_props(xr, &(top->idef), mult);
251     xr->natoms = read_first_x(oenv, &xr->traj, infile, &t, &(xr->x), xr->box);
252     xr->idef   = &(top->idef);
253     xr->oenv   = oenv;
254
255     min_max(xr);
256     calc_dihs(xr);
257
258     return top;
259 }