2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
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"
54 static const char* pp_pat[] = { "C", "N", "CA", "C", "N" };
55 #define NPP (sizeof(pp_pat) / sizeof(pp_pat[0]))
57 static bool d_comp(const t_dih& a, const t_dih& b)
59 if (a.ai[1] < b.ai[1])
63 else if (a.ai[1] == b.ai[1])
65 return a.ai[2] < b.ai[2];
74 static void calc_dihs(t_xrama* xr)
77 rvec r_ij, r_kj, r_kl, m, n;
79 gmx_rmpbc_t gpbc = nullptr;
81 gpbc = gmx_rmpbc_init(xr->idef, xr->pbcType, xr->natoms);
82 gmx_rmpbc(gpbc, xr->natoms, xr->box, xr->x);
85 for (i = 0; (i < xr->ndih); i++)
88 dd->ang = dih_angle(xr->x[dd->ai[0]],
104 gmx_bool new_data(t_xrama* xr)
106 if (!read_next_x(xr->oenv, xr->traj, &xr->t, xr->x, xr->box))
116 static int find_atom(const char* find, char*** names, int start, int nr)
120 for (i = start; (i < nr); i++)
122 if (std::strcmp(find, *names[i]) == 0)
130 static void add_xr(t_xrama* xr, const int ff[5], const t_atoms* atoms)
135 srenew(xr->dih, xr->ndih + 2LL);
136 for (i = 0; (i < 4); i++)
138 xr->dih[xr->ndih].ai[i] = ff[i];
140 for (i = 0; (i < 4); i++)
142 xr->dih[xr->ndih + 1].ai[i] = ff[i + 1];
146 srenew(xr->pp, xr->npp + 1LL);
147 xr->pp[xr->npp].iphi = xr->ndih - 2;
148 xr->pp[xr->npp].ipsi = xr->ndih - 1;
149 xr->pp[xr->npp].bShow = FALSE;
152 *atoms->resinfo[atoms->atom[ff[1]].resind].name,
153 atoms->resinfo[atoms->atom[ff[1]].resind].nr);
154 xr->pp[xr->npp].label = gmx_strdup(buf);
158 static void get_dih(t_xrama* xr, const t_atoms* atoms)
164 for (i = 0; (i < atoms->nr);)
167 for (j = 0; (j < NPP); j++)
169 if ((ff[j] = find_atom(pp_pat[j], atoms->atomname, found, atoms->nr)) == -1)
179 add_xr(xr, ff, atoms);
182 fprintf(stderr, "Found %d phi-psi combinations\n", xr->npp);
185 static void min_max(t_xrama* xr)
189 xr->amin = xr->natoms;
191 for (i = 0; (i < xr->ndih); i++)
193 for (j = 0; (j < 4); j++)
195 MSVC_DIAGNOSTIC_IGNORE(28182) // false positive in 2019 (16.5.4)
196 ai = xr->dih[i].ai[j];
197 MSVC_DIAGNOSTIC_RESET
202 else if (ai > xr->amax)
210 static void get_dih_props(t_xrama* xr, const t_idef* idef, int mult)
212 int i, ft, ftype, nra;
216 ia = idef->il[F_PDIHS].iatoms;
217 for (i = 0; (i < idef->il[F_PDIHS].nr);)
220 ftype = idef->functype[ft];
221 nra = interaction_function[ftype].nratoms;
223 if (ftype != F_PDIHS)
225 gmx_incons("ftype is not a dihedral");
230 dd = std::lower_bound(xr->dih, xr->dih + xr->ndih, key, d_comp);
231 if (dd < xr->dih + xr->ndih && !d_comp(key, *dd))
233 dd->mult = idef->iparams[ft].pdihs.mult;
234 dd->phi0 = idef->iparams[ft].pdihs.phiA;
240 /* Fill in defaults for values not in the topology */
241 for (i = 0; (i < xr->ndih); i++)
243 if (xr->dih[i].mult == 0)
246 "Dihedral around %d,%d not found in topology. Using mult=%d\n",
250 xr->dih[i].mult = mult;
251 xr->dih[i].phi0 = 180;
257 t_topology* init_rama(gmx_output_env_t* oenv, const char* infile, const char* topfile, t_xrama* xr, int mult)
262 top = read_top(topfile, &xr->pbcType);
264 /*get_dih2(xr,top->idef.functype,&(top->idef.bondeds),&(top->atoms));*/
265 get_dih(xr, &(top->atoms));
266 get_dih_props(xr, &(top->idef), mult);
267 xr->natoms = read_first_x(oenv, &xr->traj, infile, &t, &(xr->x), xr->box);
268 xr->idef = &(top->idef);