Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / gmxana / nrama.cpp
blob689f6c17ee929fda6692e54c6c39336b6bda30b4
1 /*
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,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.
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.
37 #include "gmxpre.h"
39 #include "nrama.h"
41 #include <cstdlib>
42 #include <cstring>
44 #include <algorithm>
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"
53 static const char *pp_pat[] = { "C", "N", "CA", "C", "N" };
54 #define NPP (sizeof(pp_pat)/sizeof(pp_pat[0]))
56 static bool d_comp(const t_dih &a, const t_dih &b)
58 if (a.ai[1] < b.ai[1])
60 return true;
62 else if (a.ai[1] == b.ai[1])
64 return a.ai[2] < b.ai[2];
66 else
68 return false;
73 static void calc_dihs(t_xrama *xr)
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;
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);
84 for (i = 0; (i < xr->ndih); i++)
86 dd = &(xr->dih[i]);
87 dd->ang = dih_angle(xr->x[dd->ai[0]], xr->x[dd->ai[1]],
88 xr->x[dd->ai[2]], xr->x[dd->ai[3]],
89 nullptr,
90 r_ij, r_kj, r_kl, m, n, &t1, &t2, &t3);
94 gmx_bool new_data(t_xrama *xr)
96 if (!read_next_x(xr->oenv, xr->traj, &xr->t, xr->x, xr->box))
98 return FALSE;
101 calc_dihs(xr);
103 return TRUE;
106 static int find_atom(const char *find, char ***names, int start, int nr)
108 int i;
110 for (i = start; (i < nr); i++)
112 if (std::strcmp(find, *names[i]) == 0)
114 return i;
117 return -1;
120 static void add_xr(t_xrama *xr, const int ff[5], const t_atoms *atoms)
122 char buf[12];
123 int i;
125 srenew(xr->dih, xr->ndih+2);
126 for (i = 0; (i < 4); i++)
128 xr->dih[xr->ndih].ai[i] = ff[i];
130 for (i = 0; (i < 4); i++)
132 xr->dih[xr->ndih+1].ai[i] = ff[i+1];
134 xr->ndih += 2;
136 srenew(xr->pp, xr->npp+1);
137 xr->pp[xr->npp].iphi = xr->ndih-2;
138 xr->pp[xr->npp].ipsi = xr->ndih-1;
139 xr->pp[xr->npp].bShow = FALSE;
140 sprintf(buf, "%s-%d", *atoms->resinfo[atoms->atom[ff[1]].resind].name,
141 atoms->resinfo[atoms->atom[ff[1]].resind].nr);
142 xr->pp[xr->npp].label = gmx_strdup(buf);
143 xr->npp++;
146 static void get_dih(t_xrama *xr, const t_atoms *atoms)
148 int found, ff[NPP];
149 int i;
150 size_t j;
152 for (i = 0; (i < atoms->nr); )
154 found = i;
155 for (j = 0; (j < NPP); j++)
157 if ((ff[j] = find_atom(pp_pat[j], atoms->atomname, found, atoms->nr)) == -1)
159 break;
161 found = ff[j]+1;
163 if (j != NPP)
165 break;
167 add_xr(xr, ff, atoms);
168 i = ff[0]+1;
170 fprintf(stderr, "Found %d phi-psi combinations\n", xr->npp);
173 static void min_max(t_xrama *xr)
175 int ai, i, j;
177 xr->amin = xr->natoms;
178 xr->amax = 0;
179 for (i = 0; (i < xr->ndih); i++)
181 for (j = 0; (j < 4); j++)
183 ai = xr->dih[i].ai[j];
184 if (ai < xr->amin)
186 xr->amin = ai;
188 else if (ai > xr->amax)
190 xr->amax = ai;
196 static void get_dih_props(t_xrama *xr, const t_idef *idef, int mult)
198 int i, ft, ftype, nra;
199 t_iatom *ia;
200 t_dih *dd, key;
202 ia = idef->il[F_PDIHS].iatoms;
203 for (i = 0; (i < idef->il[F_PDIHS].nr); )
205 ft = ia[0];
206 ftype = idef->functype[ft];
207 nra = interaction_function[ftype].nratoms;
209 if (ftype != F_PDIHS)
211 gmx_incons("ftype is not a dihedral");
214 key.ai[1] = ia[2];
215 key.ai[2] = ia[3];
216 dd = std::lower_bound(xr->dih, xr->dih+xr->ndih, key, d_comp);
217 if (dd < xr->dih+xr->ndih && !d_comp(key, *dd))
219 dd->mult = idef->iparams[ft].pdihs.mult;
220 dd->phi0 = idef->iparams[ft].pdihs.phiA;
223 i += nra+1;
224 ia += nra+1;
226 /* Fill in defaults for values not in the topology */
227 for (i = 0; (i < xr->ndih); i++)
229 if (xr->dih[i].mult == 0)
231 fprintf(stderr,
232 "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;
242 t_topology *init_rama(gmx_output_env_t *oenv, const char *infile,
243 const char *topfile, t_xrama *xr, int mult)
245 t_topology *top;
246 real t;
248 top = read_top(topfile, &xr->ePBC);
250 /*get_dih2(xr,top->idef.functype,&(top->idef.bondeds),&(top->atoms));*/
251 get_dih(xr, &(top->atoms));
252 get_dih_props(xr, &(top->idef), mult);
253 xr->natoms = read_first_x(oenv, &xr->traj, infile, &t, &(xr->x), xr->box);
254 xr->idef = &(top->idef);
255 xr->oenv = oenv;
257 min_max(xr);
258 calc_dihs(xr);
260 return top;