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, 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.
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.
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.
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.
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.
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.
43 #include "gromacs/listed-forces/bonded.h"
44 #include "gromacs/pbcutil/rmpbc.h"
45 #include "gromacs/utility/cstringutil.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/futil.h"
48 #include "gromacs/utility/smalloc.h"
50 static const char *pp_pat[] = { "C", "N", "CA", "C", "N" };
51 #define NPP (sizeof(pp_pat)/sizeof(pp_pat[0]))
53 static int d_comp(const void *a, const void *b)
60 if (da->ai[1] < db->ai[1])
64 else if (da->ai[1] == db->ai[1])
66 return (da->ai[2] - db->ai[2]);
75 static void calc_dihs(t_xrama *xr)
78 rvec r_ij, r_kj, r_kl, m, n;
81 gmx_rmpbc_t gpbc = NULL;
83 gpbc = gmx_rmpbc_init(xr->idef, xr->ePBC, xr->natoms);
84 gmx_rmpbc(gpbc, xr->natoms, xr->box, xr->x);
87 for (i = 0; (i < xr->ndih); i++)
90 dd->ang = dih_angle(xr->x[dd->ai[0]], xr->x[dd->ai[1]],
91 xr->x[dd->ai[2]], xr->x[dd->ai[3]],
93 r_ij, r_kj, r_kl, m, n, &sign, &t1, &t2, &t3);
97 gmx_bool new_data(t_xrama *xr)
99 if (!read_next_x(xr->oenv, xr->traj, &xr->t, xr->x, xr->box))
109 static int find_atom(const char *find, char ***names, int start, int nr)
113 for (i = start; (i < nr); i++)
115 if (strcmp(find, *names[i]) == 0)
123 static void add_xr(t_xrama *xr, int ff[5], t_atoms *atoms)
128 srenew(xr->dih, xr->ndih+2);
129 for (i = 0; (i < 4); i++)
131 xr->dih[xr->ndih].ai[i] = ff[i];
133 for (i = 0; (i < 4); i++)
135 xr->dih[xr->ndih+1].ai[i] = ff[i+1];
139 srenew(xr->pp, xr->npp+1);
140 xr->pp[xr->npp].iphi = xr->ndih-2;
141 xr->pp[xr->npp].ipsi = xr->ndih-1;
142 xr->pp[xr->npp].bShow = FALSE;
143 sprintf(buf, "%s-%d", *atoms->resinfo[atoms->atom[ff[1]].resind].name,
144 atoms->resinfo[atoms->atom[ff[1]].resind].nr);
145 xr->pp[xr->npp].label = gmx_strdup(buf);
149 static void get_dih(t_xrama *xr, t_atoms *atoms)
155 for (i = 0; (i < atoms->nr); )
158 for (j = 0; (j < NPP); j++)
160 if ((ff[j] = find_atom(pp_pat[j], atoms->atomname, found, atoms->nr)) == -1)
170 add_xr(xr, ff, atoms);
173 fprintf(stderr, "Found %d phi-psi combinations\n", xr->npp);
176 static int search_ff(int thisff[NPP], int ndih, int **ff)
179 gmx_bool bFound = FALSE;
181 for (j = 0; (j < ndih); j++)
184 for (k = 1; (k <= 3); k++)
186 bFound = bFound && (thisff[k] == ff[j][k]);
192 ff[j][4] = thisff[4];
196 ff[j][0] = thisff[0];
203 for (k = 0; (k < 5); k++)
205 ff[ndih][k] = thisff[k];
212 static void min_max(t_xrama *xr)
216 xr->amin = xr->natoms;
218 for (i = 0; (i < xr->ndih); i++)
220 for (j = 0; (j < 4); j++)
222 ai = xr->dih[i].ai[j];
227 else if (ai > xr->amax)
235 static void get_dih_props(t_xrama *xr, t_idef *idef, int mult)
237 int i, ft, ftype, nra;
241 ia = idef->il[F_PDIHS].iatoms;
242 for (i = 0; (i < idef->il[F_PDIHS].nr); )
245 ftype = idef->functype[ft];
246 nra = interaction_function[ftype].nratoms;
248 if (ftype != F_PDIHS)
250 gmx_incons("ftype is not a dihedral");
255 if ((dd = (t_dih *)bsearch(&key, xr->dih, xr->ndih, (size_t)sizeof(key), d_comp))
258 dd->mult = idef->iparams[ft].pdihs.mult;
259 dd->phi0 = idef->iparams[ft].pdihs.phiA;
265 /* Fill in defaults for values not in the topology */
266 for (i = 0; (i < xr->ndih); i++)
268 if (xr->dih[i].mult == 0)
271 "Dihedral around %d,%d not found in topology. Using mult=%d\n",
272 xr->dih[i].ai[1], xr->dih[i].ai[2], mult);
273 xr->dih[i].mult = mult;
274 xr->dih[i].phi0 = 180;
281 t_topology *init_rama(const output_env_t oenv, const char *infile,
282 const char *topfile, t_xrama *xr, int mult)
287 top = read_top(topfile, &xr->ePBC);
289 /*get_dih2(xr,top->idef.functype,&(top->idef.bondeds),&(top->atoms));*/
290 get_dih(xr, &(top->atoms));
291 get_dih_props(xr, &(top->idef), mult);
292 xr->natoms = read_first_x(oenv, &xr->traj, infile, &t, &(xr->x), xr->box);
293 xr->idef = &(top->idef);