3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
44 #include "gmx_fatal.h"
47 t_dlist *mk_dlist(FILE *log,
48 t_atoms *atoms, int *nlist,
49 gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bHChi,
50 int maxchi, int r0, gmx_residuetype_t rt)
52 int ires, i, j, k, ii;
54 int nl = 0, nc[edMax];
58 snew(dl, atoms->nres+1);
59 prev.C = prev.Cn[1] = -1; /* Keep the compiler quiet */
60 for (i = 0; (i < edMax); i++)
68 ires = atoms->atom[i].resind;
70 /* Initiate all atom numbers to -1 */
71 atm.minC = atm.H = atm.N = atm.C = atm.O = atm.minCalpha = -1;
72 for (j = 0; (j < MAXCHI+3); j++)
77 /* Look for atoms in this residue */
78 /* maybe should allow for chis to hydrogens? */
79 while ((i < atoms->nr) && (atoms->atom[i].resind == ires))
81 if ((strcmp(*(atoms->atomname[i]), "H") == 0) ||
82 (strcmp(*(atoms->atomname[i]), "H1") == 0) ||
83 (strcmp(*(atoms->atomname[i]), "HN") == 0) )
87 else if (strcmp(*(atoms->atomname[i]), "N") == 0)
91 else if (strcmp(*(atoms->atomname[i]), "C") == 0)
95 else if ((strcmp(*(atoms->atomname[i]), "O") == 0) ||
96 (strcmp(*(atoms->atomname[i]), "O1") == 0))
100 else if (strcmp(*(atoms->atomname[i]), "CA") == 0)
104 else if (strcmp(*(atoms->atomname[i]), "CB") == 0)
108 else if ((strcmp(*(atoms->atomname[i]), "CG") == 0) ||
109 (strcmp(*(atoms->atomname[i]), "CG1") == 0) ||
110 (strcmp(*(atoms->atomname[i]), "OG") == 0) ||
111 (strcmp(*(atoms->atomname[i]), "OG1") == 0) ||
112 (strcmp(*(atoms->atomname[i]), "SG") == 0))
116 else if ((strcmp(*(atoms->atomname[i]), "CD") == 0) ||
117 (strcmp(*(atoms->atomname[i]), "CD1") == 0) ||
118 (strcmp(*(atoms->atomname[i]), "SD") == 0) ||
119 (strcmp(*(atoms->atomname[i]), "OD1") == 0) ||
120 (strcmp(*(atoms->atomname[i]), "ND1") == 0))
124 /* by grs - split the Cn[4] into 2 bits to check allowing dih to H */
125 else if (bHChi && ((strcmp(*(atoms->atomname[i]), "HG") == 0) ||
126 (strcmp(*(atoms->atomname[i]), "HG1") == 0)) )
130 else if ((strcmp(*(atoms->atomname[i]), "CE") == 0) ||
131 (strcmp(*(atoms->atomname[i]), "CE1") == 0) ||
132 (strcmp(*(atoms->atomname[i]), "OE1") == 0) ||
133 (strcmp(*(atoms->atomname[i]), "NE") == 0))
137 else if ((strcmp(*(atoms->atomname[i]), "CZ") == 0) ||
138 (strcmp(*(atoms->atomname[i]), "NZ") == 0))
142 /* HChi flag here too */
143 else if (bHChi && (strcmp(*(atoms->atomname[i]), "NH1") == 0))
150 thisres = *(atoms->resinfo[ires].name);
152 /* added by grs - special case for aromatics, whose chis above 2 are
153 not real and produce rubbish output - so set back to -1 */
154 if (strcmp(thisres, "PHE") == 0 ||
155 strcmp(thisres, "TYR") == 0 ||
156 strcmp(thisres, "PTR") == 0 ||
157 strcmp(thisres, "TRP") == 0 ||
158 strcmp(thisres, "HIS") == 0 ||
159 strcmp(thisres, "HISA") == 0 ||
160 strcmp(thisres, "HISB") == 0)
162 for (ii = 5; ii <= 7; ii++)
167 /* end fixing aromatics */
169 /* Special case for Pro, has no H */
170 if (strcmp(thisres, "PRO") == 0)
174 /* Carbon from previous residue */
179 /* Alpha-carbon from previous residue */
180 if (prev.Cn[1] != -1)
182 atm.minCalpha = prev.Cn[1];
186 /* Check how many dihedrals we have */
187 if ((atm.N != -1) && (atm.Cn[1] != -1) && (atm.C != -1) &&
188 (atm.O != -1) && ((atm.H != -1) || (atm.minC != -1)))
190 dl[nl].resnr = ires+1;
192 dl[nl].atm.Cn[0] = atm.N;
193 if ((atm.Cn[3] != -1) && (atm.Cn[2] != -1) && (atm.Cn[1] != -1))
217 if ((atm.minC != -1) && (atm.minCalpha != -1))
221 dl[nl].index = gmx_residuetype_get_index(rt, thisres);
223 sprintf(dl[nl].name, "%s%d", thisres, ires+r0);
228 fprintf(debug, "Could not find N atom but could find other atoms"
229 " in residue %s%d\n", thisres, ires+r0);
232 fprintf(stderr, "\n");
234 fprintf(log, "There are %d residues with dihedrals\n", nl);
246 for (i = 0; (i < maxchi); i++)
251 fprintf(log, "There are %d dihedrals\n", j);
252 fprintf(log, "Dihedral: ");
255 fprintf(log, " Phi ");
259 fprintf(log, " Psi ");
263 for (i = 0; (i < maxchi); i++)
265 fprintf(log, "Chi%d ", i+1);
268 fprintf(log, "\nNumber: ");
271 fprintf(log, "%4d ", nl);
275 fprintf(log, "%4d ", nl);
279 for (i = 0; (i < maxchi); i++)
281 fprintf(log, "%4d ", nc[i]);
291 gmx_bool has_dihedral(int Dih, t_dlist *dl)
299 b = ((dl->atm.H != -1) && (dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1));
302 b = ((dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1) && (dl->atm.O != -1));
305 b = ((dl->atm.minCalpha != -1) && (dl->atm.minC != -1) && (dl->atm.N != -1) && (dl->atm.Cn[1] != -1));
314 b = ((dl->atm.Cn[ddd] != -1) && (dl->atm.Cn[ddd+1] != -1) &&
315 (dl->atm.Cn[ddd+2] != -1) && (dl->atm.Cn[ddd+3] != -1));
318 pr_dlist(stdout, 1, dl, 1, 0, TRUE, TRUE, TRUE, TRUE, MAXCHI);
319 gmx_fatal(FARGS, "Non existant dihedral %d in file %s, line %d",
320 Dih, __FILE__, __LINE__);
325 static void pr_one_ro(FILE *fp, t_dlist *dl, int nDih, real gmx_unused dt)
328 for (k = 0; k < NROT; k++)
330 fprintf(fp, " %6.2f", dl->rot_occ[nDih][k]);
335 static void pr_ntr_s2(FILE *fp, t_dlist *dl, int nDih, real dt)
337 fprintf(fp, " %6.2f %6.2f\n", (dt == 0) ? 0 : dl->ntr[nDih]/dt, dl->S2[nDih]);
340 void pr_dlist(FILE *fp, int nl, t_dlist dl[], real dt, int printtype,
341 gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bOmega, int maxchi)
345 void (*pr_props)(FILE *, t_dlist *, int, real);
347 /* Analysis of dihedral transitions etc */
349 if (printtype == edPrintST)
351 pr_props = pr_ntr_s2;
352 fprintf(stderr, "Now printing out transitions and OPs...\n");
356 pr_props = pr_one_ro;
357 fprintf(stderr, "Now printing out rotamer occupancies...\n");
358 fprintf(fp, "\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX\n\n");
361 /* change atom numbers from 0 based to 1 based */
362 for (i = 0; (i < nl); i++)
364 fprintf(fp, "Residue %s\n", dl[i].name);
365 if (printtype == edPrintST)
367 fprintf(fp, " Angle [ AI, AJ, AK, AL] #tr/ns S^2D \n"
368 "--------------------------------------------\n");
372 fprintf(fp, " Angle [ AI, AJ, AK, AL] rotamers 0 g(-) t g(+)\n"
373 "--------------------------------------------\n");
377 fprintf(fp, " Phi [%5d,%5d,%5d,%5d]",
378 (dl[i].atm.H == -1) ? 1+dl[i].atm.minC : 1+dl[i].atm.H,
379 1+dl[i].atm.N, 1+dl[i].atm.Cn[1], 1+dl[i].atm.C);
380 pr_props(fp, &dl[i], edPhi, dt);
384 fprintf(fp, " Psi [%5d,%5d,%5d,%5d]", 1+dl[i].atm.N, 1+dl[i].atm.Cn[1],
385 1+dl[i].atm.C, 1+dl[i].atm.O);
386 pr_props(fp, &dl[i], edPsi, dt);
388 if (bOmega && has_dihedral(edOmega, &(dl[i])))
390 fprintf(fp, " Omega [%5d,%5d,%5d,%5d]", 1+dl[i].atm.minCalpha, 1+dl[i].atm.minC,
391 1+dl[i].atm.N, 1+dl[i].atm.Cn[1]);
392 pr_props(fp, &dl[i], edOmega, dt);
394 for (Xi = 0; Xi < MAXCHI; Xi++)
396 if (bChi && (Xi < maxchi) && (dl[i].atm.Cn[Xi+3] != -1) )
398 fprintf(fp, " Chi%d[%5d,%5d,%5d,%5d]", Xi+1, 1+dl[i].atm.Cn[Xi],
399 1+dl[i].atm.Cn[Xi+1], 1+dl[i].atm.Cn[Xi+2],
400 1+dl[i].atm.Cn[Xi+3]);
401 pr_props(fp, &dl[i], Xi+edChi1, dt); /* Xi+2 was wrong here */
410 int pr_trans(FILE *fp, int nl, t_dlist dl[], real dt, int Xi)
412 /* never called at the moment */
417 fprintf(fp, "\\begin{table}[h]\n");
418 fprintf(fp, "\\caption{Number of dihedral transitions per nanosecond}\n");
419 fprintf(fp, "\\begin{tabular}{|l|l|}\n");
420 fprintf(fp, "\\hline\n");
421 fprintf(fp, "Residue\t&$\\chi_%d$\t\\\\\n", Xi+1);
422 for (i = 0; (i < nl); i++)
424 nn = dl[i].ntr[Xi]/dt;
428 fprintf(fp, "%s\t&\\HL{%d}\t\\\\\n", dl[i].name, nn);
433 fprintf(fp, "%s\t&\\%d\t\\\\\n", dl[i].name, nn);
436 fprintf(fp, "\\hline\n");
437 fprintf(fp, "\\end{tabular}\n");
438 fprintf(fp, "\\end{table}\n\n");