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.
44 #include "gromacs/listed-forces/bonded.h"
45 #include "gromacs/pbcutil/rmpbc.h"
46 #include "gromacs/utility/cstringutil.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/futil.h"
49 #include "gromacs/utility/smalloc.h"
51 static const char *pp_pat
[] = { "C", "N", "CA", "C", "N" };
52 #define NPP (sizeof(pp_pat)/sizeof(pp_pat[0]))
54 static int d_comp(const void *a
, const void *b
)
61 if (da
->ai
[1] < db
->ai
[1])
65 else if (da
->ai
[1] == db
->ai
[1])
67 return (da
->ai
[2] - db
->ai
[2]);
76 static void calc_dihs(t_xrama
*xr
)
79 rvec r_ij
, r_kj
, r_kl
, m
, n
;
82 gmx_rmpbc_t gpbc
= NULL
;
84 gpbc
= gmx_rmpbc_init(xr
->idef
, xr
->ePBC
, xr
->natoms
);
85 gmx_rmpbc(gpbc
, xr
->natoms
, xr
->box
, xr
->x
);
88 for (i
= 0; (i
< xr
->ndih
); i
++)
91 dd
->ang
= dih_angle(xr
->x
[dd
->ai
[0]], xr
->x
[dd
->ai
[1]],
92 xr
->x
[dd
->ai
[2]], xr
->x
[dd
->ai
[3]],
94 r_ij
, r_kj
, r_kl
, m
, n
, &sign
, &t1
, &t2
, &t3
);
98 gmx_bool
new_data(t_xrama
*xr
)
100 if (!read_next_x(xr
->oenv
, xr
->traj
, &xr
->t
, xr
->x
, xr
->box
))
110 static int find_atom(const char *find
, char ***names
, int start
, int nr
)
114 for (i
= start
; (i
< nr
); i
++)
116 if (std::strcmp(find
, *names
[i
]) == 0)
124 static void add_xr(t_xrama
*xr
, int ff
[5], t_atoms
*atoms
)
129 srenew(xr
->dih
, xr
->ndih
+2);
130 for (i
= 0; (i
< 4); i
++)
132 xr
->dih
[xr
->ndih
].ai
[i
] = ff
[i
];
134 for (i
= 0; (i
< 4); i
++)
136 xr
->dih
[xr
->ndih
+1].ai
[i
] = ff
[i
+1];
140 srenew(xr
->pp
, xr
->npp
+1);
141 xr
->pp
[xr
->npp
].iphi
= xr
->ndih
-2;
142 xr
->pp
[xr
->npp
].ipsi
= xr
->ndih
-1;
143 xr
->pp
[xr
->npp
].bShow
= FALSE
;
144 sprintf(buf
, "%s-%d", *atoms
->resinfo
[atoms
->atom
[ff
[1]].resind
].name
,
145 atoms
->resinfo
[atoms
->atom
[ff
[1]].resind
].nr
);
146 xr
->pp
[xr
->npp
].label
= gmx_strdup(buf
);
150 static void get_dih(t_xrama
*xr
, t_atoms
*atoms
)
156 for (i
= 0; (i
< atoms
->nr
); )
159 for (j
= 0; (j
< NPP
); j
++)
161 if ((ff
[j
] = find_atom(pp_pat
[j
], atoms
->atomname
, found
, atoms
->nr
)) == -1)
171 add_xr(xr
, ff
, atoms
);
174 fprintf(stderr
, "Found %d phi-psi combinations\n", xr
->npp
);
177 static int search_ff(int thisff
[NPP
], int ndih
, int **ff
)
180 gmx_bool bFound
= FALSE
;
182 for (j
= 0; (j
< ndih
); j
++)
185 for (k
= 1; (k
<= 3); k
++)
187 bFound
= bFound
&& (thisff
[k
] == ff
[j
][k
]);
193 ff
[j
][4] = thisff
[4];
197 ff
[j
][0] = thisff
[0];
204 for (k
= 0; (k
< 5); k
++)
206 ff
[ndih
][k
] = thisff
[k
];
213 static void min_max(t_xrama
*xr
)
217 xr
->amin
= xr
->natoms
;
219 for (i
= 0; (i
< xr
->ndih
); i
++)
221 for (j
= 0; (j
< 4); j
++)
223 ai
= xr
->dih
[i
].ai
[j
];
228 else if (ai
> xr
->amax
)
236 static void get_dih_props(t_xrama
*xr
, t_idef
*idef
, int mult
)
238 int i
, ft
, ftype
, nra
;
242 ia
= idef
->il
[F_PDIHS
].iatoms
;
243 for (i
= 0; (i
< idef
->il
[F_PDIHS
].nr
); )
246 ftype
= idef
->functype
[ft
];
247 nra
= interaction_function
[ftype
].nratoms
;
249 if (ftype
!= F_PDIHS
)
251 gmx_incons("ftype is not a dihedral");
256 if ((dd
= (t_dih
*)bsearch(&key
, xr
->dih
, xr
->ndih
, (size_t)sizeof(key
), d_comp
))
259 dd
->mult
= idef
->iparams
[ft
].pdihs
.mult
;
260 dd
->phi0
= idef
->iparams
[ft
].pdihs
.phiA
;
266 /* Fill in defaults for values not in the topology */
267 for (i
= 0; (i
< xr
->ndih
); i
++)
269 if (xr
->dih
[i
].mult
== 0)
272 "Dihedral around %d,%d not found in topology. Using mult=%d\n",
273 xr
->dih
[i
].ai
[1], xr
->dih
[i
].ai
[2], mult
);
274 xr
->dih
[i
].mult
= mult
;
275 xr
->dih
[i
].phi0
= 180;
282 t_topology
*init_rama(gmx_output_env_t
*oenv
, const char *infile
,
283 const char *topfile
, t_xrama
*xr
, int mult
)
288 top
= read_top(topfile
, &xr
->ePBC
);
290 /*get_dih2(xr,top->idef.functype,&(top->idef.bondeds),&(top->atoms));*/
291 get_dih(xr
, &(top
->atoms
));
292 get_dih_props(xr
, &(top
->idef
), mult
);
293 xr
->natoms
= read_first_x(oenv
, &xr
->traj
, infile
, &t
, &(xr
->x
), xr
->box
);
294 xr
->idef
= &(top
->idef
);