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, 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.
46 #include "gromacs/topology/residuetypes.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/smalloc.h"
50 t_dlist *mk_dlist(FILE *log,
51 t_atoms *atoms, int *nlist,
52 gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bHChi,
53 int maxchi, int r0, gmx_residuetype_t rt)
55 int ires, i, j, k, ii;
57 int nl = 0, nc[edMax];
61 snew(dl, atoms->nres+1);
62 prev.C = prev.Cn[1] = -1; /* Keep the compiler quiet */
63 for (i = 0; (i < edMax); i++)
71 ires = atoms->atom[i].resind;
73 /* Initiate all atom numbers to -1 */
74 atm.minC = atm.H = atm.N = atm.C = atm.O = atm.minCalpha = -1;
75 for (j = 0; (j < MAXCHI+3); j++)
80 /* Look for atoms in this residue */
81 /* maybe should allow for chis to hydrogens? */
82 while ((i < atoms->nr) && (atoms->atom[i].resind == ires))
84 if ((strcmp(*(atoms->atomname[i]), "H") == 0) ||
85 (strcmp(*(atoms->atomname[i]), "H1") == 0) ||
86 (strcmp(*(atoms->atomname[i]), "HN") == 0) )
90 else if (strcmp(*(atoms->atomname[i]), "N") == 0)
94 else if (strcmp(*(atoms->atomname[i]), "C") == 0)
98 else if ((strcmp(*(atoms->atomname[i]), "O") == 0) ||
99 (strcmp(*(atoms->atomname[i]), "O1") == 0) ||
100 (strcmp(*(atoms->atomname[i]), "OC1") == 0) ||
101 (strcmp(*(atoms->atomname[i]), "OT1") == 0))
105 else if (strcmp(*(atoms->atomname[i]), "CA") == 0)
109 else if (strcmp(*(atoms->atomname[i]), "CB") == 0)
113 else if ((strcmp(*(atoms->atomname[i]), "CG") == 0) ||
114 (strcmp(*(atoms->atomname[i]), "CG1") == 0) ||
115 (strcmp(*(atoms->atomname[i]), "OG") == 0) ||
116 (strcmp(*(atoms->atomname[i]), "OG1") == 0) ||
117 (strcmp(*(atoms->atomname[i]), "SG") == 0))
121 else if ((strcmp(*(atoms->atomname[i]), "CD") == 0) ||
122 (strcmp(*(atoms->atomname[i]), "CD1") == 0) ||
123 (strcmp(*(atoms->atomname[i]), "SD") == 0) ||
124 (strcmp(*(atoms->atomname[i]), "OD1") == 0) ||
125 (strcmp(*(atoms->atomname[i]), "ND1") == 0))
129 /* by grs - split the Cn[4] into 2 bits to check allowing dih to H */
130 else if (bHChi && ((strcmp(*(atoms->atomname[i]), "HG") == 0) ||
131 (strcmp(*(atoms->atomname[i]), "HG1") == 0)) )
135 else if ((strcmp(*(atoms->atomname[i]), "CE") == 0) ||
136 (strcmp(*(atoms->atomname[i]), "CE1") == 0) ||
137 (strcmp(*(atoms->atomname[i]), "OE1") == 0) ||
138 (strcmp(*(atoms->atomname[i]), "NE") == 0))
142 else if ((strcmp(*(atoms->atomname[i]), "CZ") == 0) ||
143 (strcmp(*(atoms->atomname[i]), "NZ") == 0))
147 /* HChi flag here too */
148 else if (bHChi && (strcmp(*(atoms->atomname[i]), "NH1") == 0))
155 thisres = *(atoms->resinfo[ires].name);
157 /* added by grs - special case for aromatics, whose chis above 2 are
158 not real and produce rubbish output - so set back to -1 */
159 if (strcmp(thisres, "PHE") == 0 ||
160 strcmp(thisres, "TYR") == 0 ||
161 strcmp(thisres, "PTR") == 0 ||
162 strcmp(thisres, "TRP") == 0 ||
163 strcmp(thisres, "HIS") == 0 ||
164 strcmp(thisres, "HISA") == 0 ||
165 strcmp(thisres, "HISB") == 0)
167 for (ii = 5; ii <= 7; ii++)
172 /* end fixing aromatics */
174 /* Special case for Pro, has no H */
175 if (strcmp(thisres, "PRO") == 0)
179 /* Carbon from previous residue */
184 /* Alpha-carbon from previous residue */
185 if (prev.Cn[1] != -1)
187 atm.minCalpha = prev.Cn[1];
191 /* Check how many dihedrals we have */
192 if ((atm.N != -1) && (atm.Cn[1] != -1) && (atm.C != -1) &&
193 (atm.O != -1) && ((atm.H != -1) || (atm.minC != -1)))
195 dl[nl].resnr = ires+1;
197 dl[nl].atm.Cn[0] = atm.N;
198 if ((atm.Cn[3] != -1) && (atm.Cn[2] != -1) && (atm.Cn[1] != -1))
222 if ((atm.minC != -1) && (atm.minCalpha != -1))
226 dl[nl].index = gmx_residuetype_get_index(rt, thisres);
228 sprintf(dl[nl].name, "%s%d", thisres, ires+r0);
233 fprintf(debug, "Could not find N atom but could find other atoms"
234 " in residue %s%d\n", thisres, ires+r0);
237 fprintf(stderr, "\n");
239 fprintf(log, "There are %d residues with dihedrals\n", nl);
251 for (i = 0; (i < maxchi); i++)
256 fprintf(log, "There are %d dihedrals\n", j);
257 fprintf(log, "Dihedral: ");
260 fprintf(log, " Phi ");
264 fprintf(log, " Psi ");
268 for (i = 0; (i < maxchi); i++)
270 fprintf(log, "Chi%d ", i+1);
273 fprintf(log, "\nNumber: ");
276 fprintf(log, "%4d ", nl);
280 fprintf(log, "%4d ", nl);
284 for (i = 0; (i < maxchi); i++)
286 fprintf(log, "%4d ", nc[i]);
296 gmx_bool has_dihedral(int Dih, t_dlist *dl)
304 b = ((dl->atm.H != -1) && (dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1));
307 b = ((dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1) && (dl->atm.O != -1));
310 b = ((dl->atm.minCalpha != -1) && (dl->atm.minC != -1) && (dl->atm.N != -1) && (dl->atm.Cn[1] != -1));
319 b = ((dl->atm.Cn[ddd] != -1) && (dl->atm.Cn[ddd+1] != -1) &&
320 (dl->atm.Cn[ddd+2] != -1) && (dl->atm.Cn[ddd+3] != -1));
323 pr_dlist(stdout, 1, dl, 1, 0, TRUE, TRUE, TRUE, TRUE, MAXCHI);
324 gmx_fatal(FARGS, "Non existant dihedral %d in file %s, line %d",
325 Dih, __FILE__, __LINE__);
330 static void pr_one_ro(FILE *fp, t_dlist *dl, int nDih, real gmx_unused dt)
333 for (k = 0; k < NROT; k++)
335 fprintf(fp, " %6.2f", dl->rot_occ[nDih][k]);
340 static void pr_ntr_s2(FILE *fp, t_dlist *dl, int nDih, real dt)
342 fprintf(fp, " %6.2f %6.2f\n", (dt == 0) ? 0 : dl->ntr[nDih]/dt, dl->S2[nDih]);
345 void pr_dlist(FILE *fp, int nl, t_dlist dl[], real dt, int printtype,
346 gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bOmega, int maxchi)
350 void (*pr_props)(FILE *, t_dlist *, int, real);
352 /* Analysis of dihedral transitions etc */
354 if (printtype == edPrintST)
356 pr_props = pr_ntr_s2;
357 fprintf(stderr, "Now printing out transitions and OPs...\n");
361 pr_props = pr_one_ro;
362 fprintf(stderr, "Now printing out rotamer occupancies...\n");
363 fprintf(fp, "\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX\n\n");
366 /* change atom numbers from 0 based to 1 based */
367 for (i = 0; (i < nl); i++)
369 fprintf(fp, "Residue %s\n", dl[i].name);
370 if (printtype == edPrintST)
372 fprintf(fp, " Angle [ AI, AJ, AK, AL] #tr/ns S^2D \n"
373 "--------------------------------------------\n");
377 fprintf(fp, " Angle [ AI, AJ, AK, AL] rotamers 0 g(-) t g(+)\n"
378 "--------------------------------------------\n");
382 fprintf(fp, " Phi [%5d,%5d,%5d,%5d]",
383 (dl[i].atm.H == -1) ? 1+dl[i].atm.minC : 1+dl[i].atm.H,
384 1+dl[i].atm.N, 1+dl[i].atm.Cn[1], 1+dl[i].atm.C);
385 pr_props(fp, &dl[i], edPhi, dt);
389 fprintf(fp, " Psi [%5d,%5d,%5d,%5d]", 1+dl[i].atm.N, 1+dl[i].atm.Cn[1],
390 1+dl[i].atm.C, 1+dl[i].atm.O);
391 pr_props(fp, &dl[i], edPsi, dt);
393 if (bOmega && has_dihedral(edOmega, &(dl[i])))
395 fprintf(fp, " Omega [%5d,%5d,%5d,%5d]", 1+dl[i].atm.minCalpha, 1+dl[i].atm.minC,
396 1+dl[i].atm.N, 1+dl[i].atm.Cn[1]);
397 pr_props(fp, &dl[i], edOmega, dt);
399 for (Xi = 0; Xi < MAXCHI; Xi++)
401 if (bChi && (Xi < maxchi) && (dl[i].atm.Cn[Xi+3] != -1) )
403 fprintf(fp, " Chi%d[%5d,%5d,%5d,%5d]", Xi+1, 1+dl[i].atm.Cn[Xi],
404 1+dl[i].atm.Cn[Xi+1], 1+dl[i].atm.Cn[Xi+2],
405 1+dl[i].atm.Cn[Xi+3]);
406 pr_props(fp, &dl[i], Xi+edChi1, dt); /* Xi+2 was wrong here */
415 int pr_trans(FILE *fp, int nl, t_dlist dl[], real dt, int Xi)
417 /* never called at the moment */
422 fprintf(fp, "\\begin{table}[h]\n");
423 fprintf(fp, "\\caption{Number of dihedral transitions per nanosecond}\n");
424 fprintf(fp, "\\begin{tabular}{|l|l|}\n");
425 fprintf(fp, "\\hline\n");
426 fprintf(fp, "Residue\t&$\\chi_%d$\t\\\\\n", Xi+1);
427 for (i = 0; (i < nl); i++)
429 nn = dl[i].ntr[Xi]/dt;
433 fprintf(fp, "%s\t&\\HL{%d}\t\\\\\n", dl[i].name, nn);
438 fprintf(fp, "%s\t&\\%d\t\\\\\n", dl[i].name, nn);
441 fprintf(fp, "\\hline\n");
442 fprintf(fp, "\\end{tabular}\n");
443 fprintf(fp, "\\end{table}\n\n");