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.
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/correlationfunctions/autocorr.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/matio.h"
50 #include "gromacs/fileio/pdbio.h"
51 #include "gromacs/fileio/xvgr.h"
52 #include "gromacs/gmxana/gmx_ana.h"
53 #include "gromacs/gmxana/gstat.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/utilities.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/topology/residuetypes.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/arraysize.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/smalloc.h"
65 static gmx_bool
bAllowed(real phi
, real psi
)
67 static const char *map
[] = {
68 "1100000000000000001111111000000000001111111111111111111111111",
69 "1100000000000000001111110000000000011111111111111111111111111",
70 "1100000000000000001111110000000000011111111111111111111111111",
71 "1100000000000000001111100000000000111111111111111111111111111",
72 "1100000000000000001111100000000000111111111111111111111111111",
73 "1100000000000000001111100000000001111111111111111111111111111",
74 "1100000000000000001111100000000001111111111111111111111111111",
75 "1100000000000000001111100000000011111111111111111111111111111",
76 "1110000000000000001111110000000111111111111111111111111111111",
77 "1110000000000000001111110000001111111111111111111111111111111",
78 "1110000000000000001111111000011111111111111111111111111111111",
79 "1110000000000000001111111100111111111111111111111111111111111",
80 "1110000000000000001111111111111111111111111111111111111111111",
81 "1110000000000000001111111111111111111111111111111111111111111",
82 "1110000000000000001111111111111111111111111111111111111111111",
83 "1110000000000000001111111111111111111111111111111111111111111",
84 "1110000000000000001111111111111110011111111111111111111111111",
85 "1110000000000000001111111111111100000111111111111111111111111",
86 "1110000000000000001111111111111000000000001111111111111111111",
87 "1100000000000000001111111111110000000000000011111111111111111",
88 "1100000000000000001111111111100000000000000011111111111111111",
89 "1000000000000000001111111111000000000000000001111111111111110",
90 "0000000000000000001111111110000000000000000000111111111111100",
91 "0000000000000000000000000000000000000000000000000000000000000",
92 "0000000000000000000000000000000000000000000000000000000000000",
93 "0000000000000000000000000000000000000000000000000000000000000",
94 "0000000000000000000000000000000000000000000000000000000000000",
95 "0000000000000000000000000000000000000000000000000000000000000",
96 "0000000000000000000000000000000000000000000000000000000000000",
97 "0000000000000000000000000000000000000000000000000000000000000",
98 "0000000000000000000000000000000000000000000000000000000000000",
99 "0000000000000000000000000000000000000000000000000000000000000",
100 "0000000000000000000000000000000000000000000000000000000000000",
101 "0000000000000000000000000000000000000000000000000000000000000",
102 "0000000000000000000000000000000000000000000000000000000000000",
103 "0000000000000000000000000000000000000000000000000000000000000",
104 "0000000000000000000000000000000000000000000000000000000000000",
105 "0000000000000000000000000000000000000000000000000000000000000",
106 "0000000000000000000000000000000000111111111111000000000000000",
107 "1100000000000000000000000000000001111111111111100000000000111",
108 "1100000000000000000000000000000001111111111111110000000000111",
109 "0000000000000000000000000000000000000000000000000000000000000",
110 "0000000000000000000000000000000000000000000000000000000000000",
111 "0000000000000000000000000000000000000000000000000000000000000",
112 "0000000000000000000000000000000000000000000000000000000000000",
113 "0000000000000000000000000000000000000000000000000000000000000",
114 "0000000000000000000000000000000000000000000000000000000000000",
115 "0000000000000000000000000000000000000000000000000000000000000",
116 "0000000000000000000000000000000000000000000000000000000000000",
117 "0000000000000000000000000000000000000000000000000000000000000",
118 "0000000000000000000000000000000000000000000000000000000000000",
119 "0000000000000000000000000000000000000000000000000000000000000",
120 "0000000000000000000000000000000000000000000000000000000000000",
121 "0000000000000000000000000000000000000000000000000000000000000",
122 "0000000000000000000000000000000000000000000000000000000000000",
123 "0000000000000000000000000000000000000000000000000000000000000",
124 "0000000000000000000000000000000000000000000000000000000000000",
125 "0000000000000000000000000000000000000000000000000000000000000",
126 "0000000000000000000000000000000000000000000000000000000000000",
127 "0000000000000000000000000000000000000000000000000000000000000",
128 "0000000000000000000000000000000000000000000000000000000000000"
130 #define NPP asize(map)
133 #define INDEX(ppp) (((static_cast<int> (360+ppp*RAD2DEG)) % 360)/6)
138 return (map
[x
][y
] == '1') ? TRUE
: FALSE
;
141 int *make_chi_ind(int nl
, t_dlist dl
[], int *ndih
)
146 /* There are nl residues with max edMax dihedrals with 4 atoms each */
147 snew(id
, nl
*edMax
*4);
150 for (i
= 0; (i
< nl
); i
++)
152 /* Phi, fake the first one */
153 dl
[i
].j0
[edPhi
] = n
/4;
154 if (dl
[i
].atm
.minC
>= 0)
156 id
[n
++] = dl
[i
].atm
.minC
;
160 id
[n
++] = dl
[i
].atm
.H
;
162 id
[n
++] = dl
[i
].atm
.N
;
163 id
[n
++] = dl
[i
].atm
.Cn
[1];
164 id
[n
++] = dl
[i
].atm
.C
;
166 for (i
= 0; (i
< nl
); i
++)
168 /* Psi, fake the last one */
169 dl
[i
].j0
[edPsi
] = n
/4;
170 id
[n
++] = dl
[i
].atm
.N
;
171 id
[n
++] = dl
[i
].atm
.Cn
[1];
172 id
[n
++] = dl
[i
].atm
.C
;
175 id
[n
++] = dl
[i
+1].atm
.N
;
179 id
[n
++] = dl
[i
].atm
.O
;
182 for (i
= 1; (i
< nl
); i
++)
185 if (has_dihedral(edOmega
, &(dl
[i
])))
187 dl
[i
].j0
[edOmega
] = n
/4;
188 id
[n
++] = dl
[i
].atm
.minCalpha
;
189 id
[n
++] = dl
[i
].atm
.minC
;
190 id
[n
++] = dl
[i
].atm
.N
;
191 id
[n
++] = dl
[i
].atm
.Cn
[1];
194 for (Xi
= 0; (Xi
< MAXCHI
); Xi
++)
197 for (i
= 0; (i
< nl
); i
++)
199 if (dl
[i
].atm
.Cn
[Xi
+3] != -1)
201 dl
[i
].j0
[edChi1
+Xi
] = n
/4;
202 id
[n
++] = dl
[i
].atm
.Cn
[Xi
];
203 id
[n
++] = dl
[i
].atm
.Cn
[Xi
+1];
204 id
[n
++] = dl
[i
].atm
.Cn
[Xi
+2];
205 id
[n
++] = dl
[i
].atm
.Cn
[Xi
+3];
214 int bin(real chi
, int mult
)
218 return static_cast<int>(chi
*mult
/360.0);
222 static void do_dihcorr(const char *fn
, int nf
, int ndih
, real
**dih
, real dt
,
223 int nlist
, t_dlist dlist
[], real time
[], int maxchi
,
224 gmx_bool bPhi
, gmx_bool bPsi
, gmx_bool bChi
, gmx_bool bOmega
,
225 const gmx_output_env_t
*oenv
)
227 char name1
[256], name2
[256];
230 do_autocorr(fn
, oenv
, "Dihedral Autocorrelation Function",
231 nf
, ndih
, dih
, dt
, eacCos
, FALSE
);
234 for (i
= 0; (i
< nlist
); i
++)
238 print_one(oenv
, "corrphi", dlist
[i
].name
, "Phi ACF for", "C(t)", nf
/2, time
,
243 for (i
= 0; (i
< nlist
); i
++)
247 print_one(oenv
, "corrpsi", dlist
[i
].name
, "Psi ACF for", "C(t)", nf
/2, time
,
252 for (i
= 0; (i
< nlist
); i
++)
254 if (has_dihedral(edOmega
, &dlist
[i
]))
258 print_one(oenv
, "corromega", dlist
[i
].name
, "Omega ACF for", "C(t)",
264 for (Xi
= 0; (Xi
< maxchi
); Xi
++)
266 sprintf(name1
, "corrchi%d", Xi
+1);
267 sprintf(name2
, "Chi%d ACF for", Xi
+1);
268 for (i
= 0; (i
< nlist
); i
++)
270 if (dlist
[i
].atm
.Cn
[Xi
+3] != -1)
274 print_one(oenv
, name1
, dlist
[i
].name
, name2
, "C(t)", nf
/2, time
, dih
[j
]);
280 fprintf(stderr
, "\n");
283 static void copy_dih_data(real in
[], real out
[], int nf
, gmx_bool bLEAVE
)
285 /* if bLEAVE, do nothing to data in copying to out
286 * otherwise multiply by 180/pi to convert rad to deg */
297 for (i
= 0; (i
< nf
); i
++)
303 static void dump_em_all(int nlist
, t_dlist dlist
[], int nf
, real time
[],
304 real
**dih
, int maxchi
,
305 gmx_bool bPhi
, gmx_bool bPsi
, gmx_bool bChi
, gmx_bool bOmega
, gmx_bool bRAD
,
306 const gmx_output_env_t
*oenv
)
308 char name
[256], titlestr
[256], ystr
[256];
315 std::strcpy(ystr
, "Angle (rad)");
319 std::strcpy(ystr
, "Angle (degrees)");
324 for (i
= 0; (i
< nlist
); i
++)
326 /* grs debug printf("OK i %d j %d\n", i, j) ; */
329 copy_dih_data(dih
[j
], data
, nf
, bRAD
);
330 print_one(oenv
, "phi", dlist
[i
].name
, "\\xf\\f{}", ystr
, nf
, time
, data
);
334 for (i
= 0; (i
< nlist
); i
++)
338 copy_dih_data(dih
[j
], data
, nf
, bRAD
);
339 print_one(oenv
, "psi", dlist
[i
].name
, "\\xy\\f{}", ystr
, nf
, time
, data
);
343 for (i
= 0; (i
< nlist
); i
++)
345 if (has_dihedral(edOmega
, &(dlist
[i
])))
349 copy_dih_data(dih
[j
], data
, nf
, bRAD
);
350 print_one(oenv
, "omega", dlist
[i
].name
, "\\xw\\f{}", ystr
, nf
, time
, data
);
356 for (Xi
= 0; (Xi
< maxchi
); Xi
++)
358 for (i
= 0; (i
< nlist
); i
++)
360 if (dlist
[i
].atm
.Cn
[Xi
+3] != -1)
364 sprintf(name
, "chi%d", Xi
+1);
365 sprintf(titlestr
, "\\xc\\f{}\\s%d\\N", Xi
+1);
366 copy_dih_data(dih
[j
], data
, nf
, bRAD
);
367 print_one(oenv
, name
, dlist
[i
].name
, titlestr
, ystr
, nf
, time
, data
);
373 fprintf(stderr
, "\n");
376 static void reset_one(real dih
[], int nf
, real phase
)
380 for (j
= 0; (j
< nf
); j
++)
383 while (dih
[j
] < -M_PI
)
387 while (dih
[j
] >= M_PI
)
394 static int reset_em_all(int nlist
, t_dlist dlist
[], int nf
,
395 real
**dih
, int maxchi
)
402 for (i
= 0; (i
< nlist
); i
++)
404 if (dlist
[i
].atm
.minC
== -1)
406 reset_one(dih
[j
++], nf
, M_PI
);
410 reset_one(dih
[j
++], nf
, 0);
414 for (i
= 0; (i
< nlist
-1); i
++)
416 reset_one(dih
[j
++], nf
, 0);
418 /* last Psi is faked from O */
419 reset_one(dih
[j
++], nf
, M_PI
);
422 for (i
= 0; (i
< nlist
); i
++)
424 if (has_dihedral(edOmega
, &dlist
[i
]))
426 reset_one(dih
[j
++], nf
, 0);
429 /* Chi 1 thru maxchi */
430 for (Xi
= 0; (Xi
< maxchi
); Xi
++)
432 for (i
= 0; (i
< nlist
); i
++)
434 if (dlist
[i
].atm
.Cn
[Xi
+3] != -1)
436 reset_one(dih
[j
], nf
, 0);
441 fprintf(stderr
, "j after resetting (nr. active dihedrals) = %d\n", j
);
445 static void histogramming(FILE *log
, int nbin
, gmx_residuetype_t
*rt
,
446 int nf
, int maxchi
, real
**dih
,
447 int nlist
, t_dlist dlist
[],
449 gmx_bool bPhi
, gmx_bool bPsi
, gmx_bool bOmega
, gmx_bool bChi
,
450 gmx_bool bNormalize
, gmx_bool bSSHisto
, const char *ssdump
,
451 real bfac_max
, t_atoms
*atoms
,
452 gmx_bool bDo_jc
, const char *fn
,
453 const gmx_output_env_t
*oenv
)
455 /* also gets 3J couplings and order parameters S2 */
456 t_karplus kkkphi
[] = {
457 { "J_NHa1", 6.51, -1.76, 1.6, -M_PI
/3, 0.0, 0.0 },
458 { "J_NHa2", 6.51, -1.76, 1.6, M_PI
/3, 0.0, 0.0 },
459 { "J_HaC'", 4.0, 1.1, 0.1, 0.0, 0.0, 0.0 },
460 { "J_NHCb", 4.7, -1.5, -0.2, M_PI
/3, 0.0, 0.0 },
461 { "J_Ci-1Hai", 4.5, -1.3, -1.2, 2*M_PI
/3, 0.0, 0.0 }
463 t_karplus kkkpsi
[] = {
464 { "J_HaN", -0.88, -0.61, -0.27, M_PI
/3, 0.0, 0.0 }
466 t_karplus kkkchi1
[] = {
467 { "JHaHb2", 9.5, -1.6, 1.8, -M_PI
/3, 0, 0.0 },
468 { "JHaHb3", 9.5, -1.6, 1.8, 0, 0, 0.0 }
470 #define NKKKPHI asize(kkkphi)
471 #define NKKKPSI asize(kkkpsi)
472 #define NKKKCHI asize(kkkchi1)
473 #define NJC (NKKKPHI+NKKKPSI+NKKKCHI)
475 FILE *fp
, *ssfp
[3] = {NULL
, NULL
, NULL
};
476 const char *sss
[3] = { "sheet", "helix", "coil" };
480 int ****his_aa_ss
= NULL
;
481 int ***his_aa
, *histmp
;
482 int i
, j
, k
, m
, n
, nn
, Dih
, nres
, hindex
, angle
;
483 gmx_bool bBfac
, bOccup
;
484 char hisfile
[256], hhisfile
[256], sshisfile
[256], title
[256], *ss_str
= NULL
;
486 const char *residue_name
;
489 rt_size
= gmx_residuetype_get_size(rt
);
492 fp
= gmx_ffopen(ssdump
, "r");
493 if (1 != fscanf(fp
, "%d", &nres
))
495 gmx_fatal(FARGS
, "Error reading from file %s", ssdump
);
498 snew(ss_str
, nres
+1);
499 if (1 != fscanf(fp
, "%s", ss_str
))
501 gmx_fatal(FARGS
, "Error reading from file %s", ssdump
);
505 /* Four dimensional array... Very cool */
507 for (i
= 0; (i
< 3); i
++)
509 snew(his_aa_ss
[i
], rt_size
+1);
510 for (j
= 0; (j
<= rt_size
); j
++)
512 snew(his_aa_ss
[i
][j
], edMax
);
513 for (Dih
= 0; (Dih
< edMax
); Dih
++)
515 snew(his_aa_ss
[i
][j
][Dih
], nbin
+1);
521 for (Dih
= 0; (Dih
< edMax
); Dih
++)
523 snew(his_aa
[Dih
], rt_size
+1);
524 for (i
= 0; (i
<= rt_size
); i
++)
526 snew(his_aa
[Dih
][i
], nbin
+1);
533 for (i
= 0; (i
< nlist
); i
++)
541 for (Dih
= 0; (Dih
< NONCHI
+maxchi
); Dih
++)
543 for (i
= 0; (i
< nlist
); i
++)
545 if (((Dih
< edOmega
) ) ||
546 ((Dih
== edOmega
) && (has_dihedral(edOmega
, &(dlist
[i
])))) ||
547 ((Dih
> edOmega
) && (dlist
[i
].atm
.Cn
[Dih
-NONCHI
+3] != -1)))
549 make_histo(log
, nf
, dih
[j
], nbin
, histmp
, -M_PI
, M_PI
);
553 /* Assume there is only one structure, the first.
554 * Compute index in histogram.
556 /* Check the atoms to see whether their B-factors are low enough
557 * Check atoms to see their occupancy is 1.
559 bBfac
= bOccup
= TRUE
;
560 for (nn
= 0; (nn
< 4); nn
++, n
++)
562 bBfac
= bBfac
&& (atoms
->pdbinfo
[index
[n
]].bfac
<= bfac_max
);
563 bOccup
= bOccup
&& (atoms
->pdbinfo
[index
[n
]].occup
== 1);
565 if (bOccup
&& ((bfac_max
<= 0) || ((bfac_max
> 0) && bBfac
)))
567 hindex
= static_cast<int>(((dih
[j
][0]+M_PI
)*nbin
)/(2*M_PI
));
568 range_check(hindex
, 0, nbin
);
570 /* Assign dihedral to either of the structure determined
573 switch (ss_str
[dlist
[i
].resnr
])
576 his_aa_ss
[0][dlist
[i
].index
][Dih
][hindex
]++;
579 his_aa_ss
[1][dlist
[i
].index
][Dih
][hindex
]++;
582 his_aa_ss
[2][dlist
[i
].index
][Dih
][hindex
]++;
588 fprintf(debug
, "Res. %d has imcomplete occupancy or bfacs > %g\n",
589 dlist
[i
].resnr
, bfac_max
);
600 calc_distribution_props(nbin
, histmp
, -M_PI
, NKKKPHI
, kkkphi
, &S2
);
602 for (m
= 0; (m
< NKKKPHI
); m
++)
604 Jc
[i
][m
] = kkkphi
[m
].Jc
;
605 Jcsig
[i
][m
] = kkkphi
[m
].Jcsig
;
609 calc_distribution_props(nbin
, histmp
, -M_PI
, NKKKPSI
, kkkpsi
, &S2
);
611 for (m
= 0; (m
< NKKKPSI
); m
++)
613 Jc
[i
][NKKKPHI
+m
] = kkkpsi
[m
].Jc
;
614 Jcsig
[i
][NKKKPHI
+m
] = kkkpsi
[m
].Jcsig
;
618 calc_distribution_props(nbin
, histmp
, -M_PI
, NKKKCHI
, kkkchi1
, &S2
);
619 for (m
= 0; (m
< NKKKCHI
); m
++)
621 Jc
[i
][NKKKPHI
+NKKKPSI
+m
] = kkkchi1
[m
].Jc
;
622 Jcsig
[i
][NKKKPHI
+NKKKPSI
+m
] = kkkchi1
[m
].Jcsig
;
625 default: /* covers edOmega and higher Chis than Chi1 */
626 calc_distribution_props(nbin
, histmp
, -M_PI
, 0, NULL
, &S2
);
629 dlist
[i
].S2
[Dih
] = S2
;
631 /* Sum distribution per amino acid type as well */
632 for (k
= 0; (k
< nbin
); k
++)
634 his_aa
[Dih
][dlist
[i
].index
][k
] += histmp
[k
];
639 else /* dihed not defined */
641 dlist
[i
].S2
[Dih
] = 0.0;
647 /* Print out Jcouplings */
648 fprintf(log
, "\n *** J-Couplings from simulation (plus std. dev.) ***\n\n");
649 fprintf(log
, "Residue ");
650 for (i
= 0; (i
< NKKKPHI
); i
++)
652 fprintf(log
, "%7s SD", kkkphi
[i
].name
);
654 for (i
= 0; (i
< NKKKPSI
); i
++)
656 fprintf(log
, "%7s SD", kkkpsi
[i
].name
);
658 for (i
= 0; (i
< NKKKCHI
); i
++)
660 fprintf(log
, "%7s SD", kkkchi1
[i
].name
);
663 for (i
= 0; (i
< NJC
+1); i
++)
665 fprintf(log
, "------------");
668 for (i
= 0; (i
< nlist
); i
++)
670 fprintf(log
, "%-10s", dlist
[i
].name
);
671 for (j
= 0; (j
< NJC
); j
++)
673 fprintf(log
, " %5.2f %4.2f", Jc
[i
][j
], Jcsig
[i
][j
]);
679 /* and to -jc file... */
682 fp
= xvgropen(fn
, "\\S3\\NJ-Couplings from Karplus Equation", "Residue",
685 for (i
= 0; (i
< NKKKPHI
); i
++)
687 leg
[i
] = gmx_strdup(kkkphi
[i
].name
);
689 for (i
= 0; (i
< NKKKPSI
); i
++)
691 leg
[i
+NKKKPHI
] = gmx_strdup(kkkpsi
[i
].name
);
693 for (i
= 0; (i
< NKKKCHI
); i
++)
695 leg
[i
+NKKKPHI
+NKKKPSI
] = gmx_strdup(kkkchi1
[i
].name
);
697 xvgr_legend(fp
, NJC
, (const char**)leg
, oenv
);
698 fprintf(fp
, "%5s ", "#Res.");
699 for (i
= 0; (i
< NJC
); i
++)
701 fprintf(fp
, "%10s ", leg
[i
]);
704 for (i
= 0; (i
< nlist
); i
++)
706 fprintf(fp
, "%5d ", dlist
[i
].resnr
);
707 for (j
= 0; (j
< NJC
); j
++)
709 fprintf(fp
, " %8.3f", Jc
[i
][j
]);
714 for (i
= 0; (i
< NJC
); i
++)
719 /* finished -jc stuff */
721 snew(normhisto
, nbin
);
722 for (i
= 0; (i
< rt_size
); i
++)
724 for (Dih
= 0; (Dih
< edMax
); Dih
++)
726 /* First check whether something is in there */
727 for (j
= 0; (j
< nbin
); j
++)
729 if (his_aa
[Dih
][i
][j
] != 0)
735 ((bPhi
&& (Dih
== edPhi
)) ||
736 (bPsi
&& (Dih
== edPsi
)) ||
737 (bOmega
&& (Dih
== edOmega
)) ||
738 (bChi
&& (Dih
>= edChi1
))))
742 normalize_histo(nbin
, his_aa
[Dih
][i
], (360.0/nbin
), normhisto
);
745 residue_name
= gmx_residuetype_get_name(rt
, i
);
749 sprintf(hisfile
, "histo-phi%s", residue_name
);
750 sprintf(title
, "\\xf\\f{} Distribution for %s", residue_name
);
753 sprintf(hisfile
, "histo-psi%s", residue_name
);
754 sprintf(title
, "\\xy\\f{} Distribution for %s", residue_name
);
757 sprintf(hisfile
, "histo-omega%s", residue_name
);
758 sprintf(title
, "\\xw\\f{} Distribution for %s", residue_name
);
761 sprintf(hisfile
, "histo-chi%d%s", Dih
-NONCHI
+1, residue_name
);
762 sprintf(title
, "\\xc\\f{}\\s%d\\N Distribution for %s",
763 Dih
-NONCHI
+1, residue_name
);
765 std::strcpy(hhisfile
, hisfile
);
766 std::strcat(hhisfile
, ".xvg");
767 fp
= xvgropen(hhisfile
, title
, "Degrees", "", oenv
);
768 if (output_env_get_print_xvgr_codes(oenv
))
770 fprintf(fp
, "@ with g0\n");
772 xvgr_world(fp
, -180, 0, 180, 0.1, oenv
);
773 if (output_env_get_print_xvgr_codes(oenv
))
775 fprintf(fp
, "# this effort to set graph size fails unless you run with -autoscale none or -autoscale y flags\n");
776 fprintf(fp
, "@ xaxis tick on\n");
777 fprintf(fp
, "@ xaxis tick major 90\n");
778 fprintf(fp
, "@ xaxis tick minor 30\n");
779 fprintf(fp
, "@ xaxis ticklabel prec 0\n");
780 fprintf(fp
, "@ yaxis tick off\n");
781 fprintf(fp
, "@ yaxis ticklabel off\n");
782 fprintf(fp
, "@ type xy\n");
786 for (k
= 0; (k
< 3); k
++)
788 sprintf(sshisfile
, "%s-%s.xvg", hisfile
, sss
[k
]);
789 ssfp
[k
] = gmx_ffopen(sshisfile
, "w");
792 for (j
= 0; (j
< nbin
); j
++)
794 angle
= -180 + (360/nbin
)*j
;
797 fprintf(fp
, "%5d %10g\n", angle
, normhisto
[j
]);
801 fprintf(fp
, "%5d %10d\n", angle
, his_aa
[Dih
][i
][j
]);
805 for (k
= 0; (k
< 3); k
++)
807 fprintf(ssfp
[k
], "%5d %10d\n", angle
,
808 his_aa_ss
[k
][i
][Dih
][j
]);
812 fprintf(fp
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
816 for (k
= 0; (k
< 3); k
++)
818 fprintf(ssfp
[k
], "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
819 gmx_ffclose(ssfp
[k
]);
829 /* Four dimensional array... Very cool */
830 for (i
= 0; (i
< 3); i
++)
832 for (j
= 0; (j
<= rt_size
); j
++)
834 for (Dih
= 0; (Dih
< edMax
); Dih
++)
836 sfree(his_aa_ss
[i
][j
][Dih
]);
838 sfree(his_aa_ss
[i
][j
]);
847 static FILE *rama_file(const char *fn
, const char *title
, const char *xaxis
,
848 const char *yaxis
, const gmx_output_env_t
*oenv
)
852 fp
= xvgropen(fn
, title
, xaxis
, yaxis
, oenv
);
853 if (output_env_get_print_xvgr_codes(oenv
))
855 fprintf(fp
, "@ with g0\n");
857 xvgr_world(fp
, -180, -180, 180, 180, oenv
);
858 if (output_env_get_print_xvgr_codes(oenv
))
860 fprintf(fp
, "@ xaxis tick on\n");
861 fprintf(fp
, "@ xaxis tick major 90\n");
862 fprintf(fp
, "@ xaxis tick minor 30\n");
863 fprintf(fp
, "@ xaxis ticklabel prec 0\n");
864 fprintf(fp
, "@ yaxis tick on\n");
865 fprintf(fp
, "@ yaxis tick major 90\n");
866 fprintf(fp
, "@ yaxis tick minor 30\n");
867 fprintf(fp
, "@ yaxis ticklabel prec 0\n");
868 fprintf(fp
, "@ s0 type xy\n");
869 fprintf(fp
, "@ s0 symbol 2\n");
870 fprintf(fp
, "@ s0 symbol size 0.410000\n");
871 fprintf(fp
, "@ s0 symbol fill 1\n");
872 fprintf(fp
, "@ s0 symbol color 1\n");
873 fprintf(fp
, "@ s0 symbol linewidth 1\n");
874 fprintf(fp
, "@ s0 symbol linestyle 1\n");
875 fprintf(fp
, "@ s0 symbol center false\n");
876 fprintf(fp
, "@ s0 symbol char 0\n");
877 fprintf(fp
, "@ s0 skip 0\n");
878 fprintf(fp
, "@ s0 linestyle 0\n");
879 fprintf(fp
, "@ s0 linewidth 1\n");
880 fprintf(fp
, "@ type xy\n");
885 static void do_rama(int nf
, int nlist
, t_dlist dlist
[], real
**dih
,
886 gmx_bool bViol
, gmx_bool bRamOmega
, const gmx_output_env_t
*oenv
)
888 FILE *fp
, *gp
= NULL
;
891 int i
, j
, k
, Xi1
, Xi2
, Phi
, Psi
, Om
= 0, nlevels
;
893 real
**mat
= NULL
, phi
, psi
, omega
, axis
[NMAT
], lo
, hi
;
894 t_rgb rlo
= { 1.0, 0.0, 0.0 };
895 t_rgb rmid
= { 1.0, 1.0, 1.0 };
896 t_rgb rhi
= { 0.0, 0.0, 1.0 };
898 for (i
= 0; (i
< nlist
); i
++)
900 if ((has_dihedral(edPhi
, &(dlist
[i
]))) &&
901 (has_dihedral(edPsi
, &(dlist
[i
]))))
903 sprintf(fn
, "ramaPhiPsi%s.xvg", dlist
[i
].name
);
904 fp
= rama_file(fn
, "Ramachandran Plot",
905 "\\8f\\4 (deg)", "\\8y\\4 (deg)", oenv
);
906 bOm
= bRamOmega
&& has_dihedral(edOmega
, &(dlist
[i
]));
909 Om
= dlist
[i
].j0
[edOmega
];
911 for (j
= 0; (j
< NMAT
); j
++)
914 axis
[j
] = -180+(360*j
)/NMAT
;
919 sprintf(fn
, "violPhiPsi%s.xvg", dlist
[i
].name
);
920 gp
= gmx_ffopen(fn
, "w");
922 Phi
= dlist
[i
].j0
[edPhi
];
923 Psi
= dlist
[i
].j0
[edPsi
];
924 for (j
= 0; (j
< nf
); j
++)
926 phi
= RAD2DEG
*dih
[Phi
][j
];
927 psi
= RAD2DEG
*dih
[Psi
][j
];
928 fprintf(fp
, "%10g %10g\n", phi
, psi
);
931 fprintf(gp
, "%d\n", (bAllowed(dih
[Phi
][j
], RAD2DEG
*dih
[Psi
][j
]) == FALSE
) );
935 omega
= RAD2DEG
*dih
[Om
][j
];
936 mat
[static_cast<int>(((phi
*NMAT
)/360)+NMAT
/2)][static_cast<int>(((psi
*NMAT
)/360)+NMAT
/2)]
947 sprintf(fn
, "ramomega%s.xpm", dlist
[i
].name
);
948 fp
= gmx_ffopen(fn
, "w");
950 for (j
= 0; (j
< NMAT
); j
++)
952 for (k
= 0; (k
< NMAT
); k
++)
955 lo
= std::min(mat
[j
][k
], lo
);
956 hi
= std::max(mat
[j
][k
], hi
);
960 if (std::abs(lo
) > std::abs(hi
))
969 for (j
= 0; (j
< NMAT
); j
++)
971 for (k
= 0; (k
< NMAT
); k
++)
979 write_xpm3(fp
, 0, "Omega/Ramachandran Plot", "Deg", "Phi", "Psi",
980 NMAT
, NMAT
, axis
, axis
, mat
, lo
, 180.0, hi
, rlo
, rmid
, rhi
, &nlevels
);
982 for (j
= 0; (j
< NMAT
); j
++)
989 if ((has_dihedral(edChi1
, &(dlist
[i
]))) &&
990 (has_dihedral(edChi2
, &(dlist
[i
]))))
992 sprintf(fn
, "ramaX1X2%s.xvg", dlist
[i
].name
);
993 fp
= rama_file(fn
, "\\8c\\4\\s1\\N-\\8c\\4\\s2\\N Ramachandran Plot",
994 "\\8c\\4\\s1\\N (deg)", "\\8c\\4\\s2\\N (deg)", oenv
);
995 Xi1
= dlist
[i
].j0
[edChi1
];
996 Xi2
= dlist
[i
].j0
[edChi2
];
997 for (j
= 0; (j
< nf
); j
++)
999 fprintf(fp
, "%10g %10g\n", RAD2DEG
*dih
[Xi1
][j
], RAD2DEG
*dih
[Xi2
][j
]);
1005 fprintf(stderr
, "No chi1 & chi2 angle for %s\n", dlist
[i
].name
);
1011 static void print_transitions(const char *fn
, int maxchi
, int nlist
,
1012 t_dlist dlist
[], real dt
,
1013 const gmx_output_env_t
*oenv
)
1015 /* based on order_params below */
1019 /* must correspond with enum in pp2shift.h:38 */
1021 #define NLEG asize(leg)
1023 leg
[0] = gmx_strdup("Phi");
1024 leg
[1] = gmx_strdup("Psi");
1025 leg
[2] = gmx_strdup("Omega");
1026 leg
[3] = gmx_strdup("Chi1");
1027 leg
[4] = gmx_strdup("Chi2");
1028 leg
[5] = gmx_strdup("Chi3");
1029 leg
[6] = gmx_strdup("Chi4");
1030 leg
[7] = gmx_strdup("Chi5");
1031 leg
[8] = gmx_strdup("Chi6");
1033 /* Print order parameters */
1034 fp
= xvgropen(fn
, "Dihedral Rotamer Transitions", "Residue", "Transitions/ns",
1036 xvgr_legend(fp
, NONCHI
+maxchi
, (const char**)leg
, oenv
);
1038 fprintf(fp
, "%5s ", "#Res.");
1039 fprintf(fp
, "%10s %10s %10s ", leg
[edPhi
], leg
[edPsi
], leg
[edOmega
]);
1040 for (Xi
= 0; Xi
< maxchi
; Xi
++)
1042 fprintf(fp
, "%10s ", leg
[NONCHI
+Xi
]);
1046 for (i
= 0; (i
< nlist
); i
++)
1048 fprintf(fp
, "%5d ", dlist
[i
].resnr
);
1049 for (Dih
= 0; (Dih
< NONCHI
+maxchi
); Dih
++)
1051 fprintf(fp
, "%10.3f ", dlist
[i
].ntr
[Dih
]/dt
);
1053 /* fprintf(fp,"%12s\n",dlist[i].name); this confuses xmgrace */
1059 static void order_params(FILE *log
,
1060 const char *fn
, int maxchi
, int nlist
, t_dlist dlist
[],
1061 const char *pdbfn
, real bfac_init
,
1062 t_atoms
*atoms
, rvec x
[], int ePBC
, matrix box
,
1063 gmx_bool bPhi
, gmx_bool bPsi
, gmx_bool bChi
, const gmx_output_env_t
*oenv
)
1070 /* except for S2Min/Max, must correspond with enum in pp2shift.h:38 */
1071 const char *const_leg
[2+edMax
] = {
1072 "S2Min", "S2Max", "Phi", "Psi", "Omega",
1073 "Chi1", "Chi2", "Chi3", "Chi4", "Chi5",
1076 #define NLEG asize(leg)
1080 for (i
= 0; i
< NLEG
; i
++)
1082 leg
[i
] = gmx_strdup(const_leg
[i
]);
1085 /* Print order parameters */
1086 fp
= xvgropen(fn
, "Dihedral Order Parameters", "Residue", "S2", oenv
);
1087 xvgr_legend(fp
, 2+NONCHI
+maxchi
, const_leg
, oenv
);
1089 for (Dih
= 0; (Dih
< edMax
); Dih
++)
1094 fprintf(fp
, "%5s ", "#Res.");
1095 fprintf(fp
, "%10s %10s ", leg
[0], leg
[1]);
1096 fprintf(fp
, "%10s %10s %10s ", leg
[2+edPhi
], leg
[2+edPsi
], leg
[2+edOmega
]);
1097 for (Xi
= 0; Xi
< maxchi
; Xi
++)
1099 fprintf(fp
, "%10s ", leg
[2+NONCHI
+Xi
]);
1103 for (i
= 0; (i
< nlist
); i
++)
1107 for (Dih
= 0; (Dih
< NONCHI
+maxchi
); Dih
++)
1109 if (dlist
[i
].S2
[Dih
] != 0)
1111 if (dlist
[i
].S2
[Dih
] > S2Max
)
1113 S2Max
= dlist
[i
].S2
[Dih
];
1115 if (dlist
[i
].S2
[Dih
] < S2Min
)
1117 S2Min
= dlist
[i
].S2
[Dih
];
1120 if (dlist
[i
].S2
[Dih
] > 0.8)
1125 fprintf(fp
, "%5d ", dlist
[i
].resnr
);
1126 fprintf(fp
, "%10.3f %10.3f ", S2Min
, S2Max
);
1127 for (Dih
= 0; (Dih
< NONCHI
+maxchi
); Dih
++)
1129 fprintf(fp
, "%10.3f ", dlist
[i
].S2
[Dih
]);
1132 /* fprintf(fp,"%12s\n",dlist[i].name); this confuses xmgrace */
1140 if (NULL
== atoms
->pdbinfo
)
1142 snew(atoms
->pdbinfo
, atoms
->nr
);
1144 for (i
= 0; (i
< atoms
->nr
); i
++)
1146 atoms
->pdbinfo
[i
].bfac
= bfac_init
;
1149 for (i
= 0; (i
< nlist
); i
++)
1151 atoms
->pdbinfo
[dlist
[i
].atm
.N
].bfac
= -dlist
[i
].S2
[0]; /* Phi */
1152 atoms
->pdbinfo
[dlist
[i
].atm
.H
].bfac
= -dlist
[i
].S2
[0]; /* Phi */
1153 atoms
->pdbinfo
[dlist
[i
].atm
.C
].bfac
= -dlist
[i
].S2
[1]; /* Psi */
1154 atoms
->pdbinfo
[dlist
[i
].atm
.O
].bfac
= -dlist
[i
].S2
[1]; /* Psi */
1155 for (Xi
= 0; (Xi
< maxchi
); Xi
++) /* Chi's */
1157 if (dlist
[i
].atm
.Cn
[Xi
+3] != -1)
1159 atoms
->pdbinfo
[dlist
[i
].atm
.Cn
[Xi
+1]].bfac
= -dlist
[i
].S2
[NONCHI
+Xi
];
1164 fp
= gmx_ffopen(pdbfn
, "w");
1165 fprintf(fp
, "REMARK generated by g_chi\n");
1166 fprintf(fp
, "REMARK "
1167 "B-factor field contains negative of dihedral order parameters\n");
1168 write_pdbfile(fp
, NULL
, atoms
, x
, ePBC
, box
, ' ', 0, NULL
, TRUE
);
1169 x0
= y0
= z0
= 1000.0;
1170 for (i
= 0; (i
< atoms
->nr
); i
++)
1172 x0
= std::min(x0
, x
[i
][XX
]);
1173 y0
= std::min(y0
, x
[i
][YY
]);
1174 z0
= std::min(z0
, x
[i
][ZZ
]);
1176 x0
*= 10.0; /* nm -> angstrom */
1177 y0
*= 10.0; /* nm -> angstrom */
1178 z0
*= 10.0; /* nm -> angstrom */
1179 for (i
= 0; (i
< 10); i
++)
1181 gmx_fprintf_pdb_atomline(fp
, epdbATOM
, atoms
->nr
+1+i
, "CA", ' ', "LEG", ' ', atoms
->nres
+1, ' ',
1182 x0
, y0
, z0
+(1.2*i
), 0.0, -0.1*i
, "");
1187 fprintf(log
, "Dihedrals with S2 > 0.8\n");
1188 fprintf(log
, "Dihedral: ");
1191 fprintf(log
, " Phi ");
1195 fprintf(log
, " Psi ");
1199 for (Xi
= 0; (Xi
< maxchi
); Xi
++)
1201 fprintf(log
, " %s ", leg
[2+NONCHI
+Xi
]);
1204 fprintf(log
, "\nNumber: ");
1207 fprintf(log
, "%4d ", nh
[0]);
1211 fprintf(log
, "%4d ", nh
[1]);
1215 for (Xi
= 0; (Xi
< maxchi
); Xi
++)
1217 fprintf(log
, "%4d ", nh
[NONCHI
+Xi
]);
1222 for (i
= 0; i
< NLEG
; i
++)
1229 int gmx_chi(int argc
, char *argv
[])
1231 const char *desc
[] = {
1232 "[THISMODULE] computes [GRK]phi[grk], [GRK]psi[grk], [GRK]omega[grk],",
1233 "and [GRK]chi[grk] dihedrals for all your",
1234 "amino acid backbone and sidechains.",
1235 "It can compute dihedral angle as a function of time, and as",
1236 "histogram distributions.",
1237 "The distributions [TT](histo-(dihedral)(RESIDUE).xvg[tt]) are cumulative over all residues of each type.[PAR]",
1238 "If option [TT]-corr[tt] is given, the program will",
1239 "calculate dihedral autocorrelation functions. The function used",
1240 "is C(t) = [CHEVRON][COS][GRK]chi[grk]([GRK]tau[grk])[cos] [COS][GRK]chi[grk]([GRK]tau[grk]+t)[cos][chevron]. The use of cosines",
1241 "rather than angles themselves, resolves the problem of periodicity.",
1242 "(Van der Spoel & Berendsen (1997), Biophys. J. 72, 2032-2041).",
1243 "Separate files for each dihedral of each residue",
1244 "[TT](corr(dihedral)(RESIDUE)(nresnr).xvg[tt]) are output, as well as a",
1245 "file containing the information for all residues (argument of [TT]-corr[tt]).[PAR]",
1246 "With option [TT]-all[tt], the angles themselves as a function of time for",
1247 "each residue are printed to separate files [TT](dihedral)(RESIDUE)(nresnr).xvg[tt].",
1248 "These can be in radians or degrees.[PAR]",
1249 "A log file (argument [TT]-g[tt]) is also written. This contains",
1251 " * information about the number of residues of each type.",
1252 " * The NMR ^3J coupling constants from the Karplus equation.",
1253 " * a table for each residue of the number of transitions between ",
1254 " rotamers per nanosecond, and the order parameter S^2 of each dihedral.",
1255 " * a table for each residue of the rotamer occupancy.",
1257 "All rotamers are taken as 3-fold, except for [GRK]omega[grk] and [GRK]chi[grk] dihedrals",
1258 "to planar groups (i.e. [GRK]chi[grk][SUB]2[sub] of aromatics, Asp and Asn; [GRK]chi[grk][SUB]3[sub] of Glu",
1259 "and Gln; and [GRK]chi[grk][SUB]4[sub] of Arg), which are 2-fold. \"rotamer 0\" means ",
1260 "that the dihedral was not in the core region of each rotamer. ",
1261 "The width of the core region can be set with [TT]-core_rotamer[tt][PAR]",
1263 "The S^2 order parameters are also output to an [REF].xvg[ref] file",
1264 "(argument [TT]-o[tt] ) and optionally as a [REF].pdb[ref] file with",
1265 "the S^2 values as B-factor (argument [TT]-p[tt]). ",
1266 "The total number of rotamer transitions per timestep",
1267 "(argument [TT]-ot[tt]), the number of transitions per rotamer",
1268 "(argument [TT]-rt[tt]), and the ^3J couplings (argument [TT]-jc[tt]), ",
1269 "can also be written to [REF].xvg[ref] files. Note that the analysis",
1270 "of rotamer transitions assumes that the supplied trajectory frames",
1271 "are equally spaced in time.[PAR]",
1273 "If [TT]-chi_prod[tt] is set (and [TT]-maxchi[tt] > 0), cumulative rotamers, e.g.",
1274 "1+9([GRK]chi[grk][SUB]1[sub]-1)+3([GRK]chi[grk][SUB]2[sub]-1)+([GRK]chi[grk][SUB]3[sub]-1) (if the residue has three 3-fold ",
1275 "dihedrals and [TT]-maxchi[tt] >= 3)",
1276 "are calculated. As before, if any dihedral is not in the core region,",
1277 "the rotamer is taken to be 0. The occupancies of these cumulative ",
1278 "rotamers (starting with rotamer 0) are written to the file",
1279 "that is the argument of [TT]-cp[tt], and if the [TT]-all[tt] flag",
1280 "is given, the rotamers as functions of time",
1281 "are written to [TT]chiproduct(RESIDUE)(nresnr).xvg[tt] ",
1282 "and their occupancies to [TT]histo-chiproduct(RESIDUE)(nresnr).xvg[tt].[PAR]",
1284 "The option [TT]-r[tt] generates a contour plot of the average [GRK]omega[grk] angle",
1285 "as a function of the [GRK]phi[grk] and [GRK]psi[grk] angles, that is, in a Ramachandran plot",
1286 "the average [GRK]omega[grk] angle is plotted using color coding.",
1290 const char *bugs
[] = {
1291 "Produces MANY output files (up to about 4 times the number of residues in the protein, twice that if autocorrelation functions are calculated). Typically several hundred files are output.",
1292 "[GRK]phi[grk] and [GRK]psi[grk] dihedrals are calculated in a "
1293 "non-standard way, using H-N-CA-C for [GRK]phi[grk] instead of "
1294 "C(-)-N-CA-C, and N-CA-C-O for [GRK]psi[grk] instead of N-CA-C-N(+). "
1295 "This causes (usually small) discrepancies with the output of other "
1296 "tools like [gmx-rama].",
1297 "[TT]-r0[tt] option does not work properly",
1298 "Rotamers with multiplicity 2 are printed in [TT]chi.log[tt] as if they had multiplicity 3, with the 3rd (g(+)) always having probability 0"
1302 static int r0
= 1, ndeg
= 1, maxchi
= 2;
1303 static gmx_bool bAll
= FALSE
;
1304 static gmx_bool bPhi
= FALSE
, bPsi
= FALSE
, bOmega
= FALSE
;
1305 static real bfac_init
= -1.0, bfac_max
= 0;
1306 static const char *maxchistr
[] = { NULL
, "0", "1", "2", "3", "4", "5", "6", NULL
};
1307 static gmx_bool bRama
= FALSE
, bShift
= FALSE
, bViol
= FALSE
, bRamOmega
= FALSE
;
1308 static gmx_bool bNormHisto
= TRUE
, bChiProduct
= FALSE
, bHChi
= FALSE
, bRAD
= FALSE
, bPBC
= TRUE
;
1309 static real core_frac
= 0.5;
1311 { "-r0", FALSE
, etINT
, {&r0
},
1312 "starting residue" },
1313 { "-phi", FALSE
, etBOOL
, {&bPhi
},
1314 "Output for [GRK]phi[grk] dihedral angles" },
1315 { "-psi", FALSE
, etBOOL
, {&bPsi
},
1316 "Output for [GRK]psi[grk] dihedral angles" },
1317 { "-omega", FALSE
, etBOOL
, {&bOmega
},
1318 "Output for [GRK]omega[grk] dihedrals (peptide bonds)" },
1319 { "-rama", FALSE
, etBOOL
, {&bRama
},
1320 "Generate [GRK]phi[grk]/[GRK]psi[grk] and [GRK]chi[grk][SUB]1[sub]/[GRK]chi[grk][SUB]2[sub] Ramachandran plots" },
1321 { "-viol", FALSE
, etBOOL
, {&bViol
},
1322 "Write a file that gives 0 or 1 for violated Ramachandran angles" },
1323 { "-periodic", FALSE
, etBOOL
, {&bPBC
},
1324 "Print dihedral angles modulo 360 degrees" },
1325 { "-all", FALSE
, etBOOL
, {&bAll
},
1326 "Output separate files for every dihedral." },
1327 { "-rad", FALSE
, etBOOL
, {&bRAD
},
1328 "in angle vs time files, use radians rather than degrees."},
1329 { "-shift", FALSE
, etBOOL
, {&bShift
},
1330 "Compute chemical shifts from [GRK]phi[grk]/[GRK]psi[grk] angles" },
1331 { "-binwidth", FALSE
, etINT
, {&ndeg
},
1332 "bin width for histograms (degrees)" },
1333 { "-core_rotamer", FALSE
, etREAL
, {&core_frac
},
1334 "only the central [TT]-core_rotamer[tt]\\*(360/multiplicity) belongs to each rotamer (the rest is assigned to rotamer 0)" },
1335 { "-maxchi", FALSE
, etENUM
, {maxchistr
},
1336 "calculate first ndih [GRK]chi[grk] dihedrals" },
1337 { "-normhisto", FALSE
, etBOOL
, {&bNormHisto
},
1338 "Normalize histograms" },
1339 { "-ramomega", FALSE
, etBOOL
, {&bRamOmega
},
1340 "compute average omega as a function of [GRK]phi[grk]/[GRK]psi[grk] and plot it in an [REF].xpm[ref] plot" },
1341 { "-bfact", FALSE
, etREAL
, {&bfac_init
},
1342 "B-factor value for [REF].pdb[ref] file for atoms with no calculated dihedral order parameter"},
1343 { "-chi_prod", FALSE
, etBOOL
, {&bChiProduct
},
1344 "compute a single cumulative rotamer for each residue"},
1345 { "-HChi", FALSE
, etBOOL
, {&bHChi
},
1346 "Include dihedrals to sidechain hydrogens"},
1347 { "-bmax", FALSE
, etREAL
, {&bfac_max
},
1348 "Maximum B-factor on any of the atoms that make up a dihedral, for the dihedral angle to be considere in the statistics. Applies to database work where a number of X-Ray structures is analyzed. [TT]-bmax[tt] <= 0 means no limit." }
1352 int nlist
, idum
, nbin
;
1358 gmx_bool bChi
, bCorr
, bSSHisto
;
1359 gmx_bool bDo_rt
, bDo_oh
, bDo_ot
, bDo_jc
;
1360 real dt
= 0, traj_t_ns
;
1361 gmx_output_env_t
*oenv
;
1362 gmx_residuetype_t
*rt
;
1365 int ndih
, nactdih
, nf
;
1366 real
**dih
, *trans_frac
, *aver_angle
, *time
;
1367 int i
, **chi_lookup
, *multiplicity
;
1370 { efSTX
, "-s", NULL
, ffREAD
},
1371 { efTRX
, "-f", NULL
, ffREAD
},
1372 { efXVG
, "-o", "order", ffWRITE
},
1373 { efPDB
, "-p", "order", ffOPTWR
},
1374 { efDAT
, "-ss", "ssdump", ffOPTRD
},
1375 { efXVG
, "-jc", "Jcoupling", ffWRITE
},
1376 { efXVG
, "-corr", "dihcorr", ffOPTWR
},
1377 { efLOG
, "-g", "chi", ffWRITE
},
1378 /* add two more arguments copying from g_angle */
1379 { efXVG
, "-ot", "dihtrans", ffOPTWR
},
1380 { efXVG
, "-oh", "trhisto", ffOPTWR
},
1381 { efXVG
, "-rt", "restrans", ffOPTWR
},
1382 { efXVG
, "-cp", "chiprodhisto", ffOPTWR
}
1384 #define NFILE asize(fnm)
1389 ppa
= add_acf_pargs(&npargs
, pa
);
1390 if (!parse_common_args(&argc
, argv
, PCA_CAN_VIEW
| PCA_CAN_TIME
,
1391 NFILE
, fnm
, npargs
, ppa
, asize(desc
), desc
, asize(bugs
), bugs
,
1397 /* Handle result from enumerated type */
1398 sscanf(maxchistr
[0], "%d", &maxchi
);
1399 bChi
= (maxchi
> 0);
1401 log
= gmx_ffopen(ftp2fn(efLOG
, NFILE
, fnm
), "w");
1410 /* set some options */
1411 bDo_rt
= (opt2bSet("-rt", NFILE
, fnm
));
1412 bDo_oh
= (opt2bSet("-oh", NFILE
, fnm
));
1413 bDo_ot
= (opt2bSet("-ot", NFILE
, fnm
));
1414 bDo_jc
= (opt2bSet("-jc", NFILE
, fnm
));
1415 bCorr
= (opt2bSet("-corr", NFILE
, fnm
));
1418 fprintf(stderr
, "Will calculate autocorrelation\n");
1421 if (core_frac
> 1.0)
1423 fprintf(stderr
, "core_rotamer fraction > 1.0 ; will use 1.0\n");
1426 if (core_frac
< 0.0)
1428 fprintf(stderr
, "core_rotamer fraction < 0.0 ; will use 0.0\n");
1432 if (maxchi
> MAXCHI
)
1435 "Will only calculate first %d Chi dihedrals in stead of %d.\n",
1439 bSSHisto
= ftp2bSet(efDAT
, NFILE
, fnm
);
1442 /* Find the chi angles using atoms struct and a list of amino acids */
1445 read_tps_conf(ftp2fn(efSTX
, NFILE
, fnm
), top
, &ePBC
, &x
, NULL
, box
, FALSE
);
1446 t_atoms
&atoms
= top
->atoms
;
1447 if (atoms
.pdbinfo
== NULL
)
1449 snew(atoms
.pdbinfo
, atoms
.nr
);
1451 fprintf(log
, "Title: %s\n", *top
->name
);
1453 gmx_residuetype_init(&rt
);
1454 dlist
= mk_dlist(log
, &atoms
, &nlist
, bPhi
, bPsi
, bChi
, bHChi
, maxchi
, r0
, rt
);
1455 fprintf(stderr
, "%d residues with dihedrals found\n", nlist
);
1459 gmx_fatal(FARGS
, "No dihedrals in your structure!\n");
1462 /* Make a linear index for reading all. */
1463 index
= make_chi_ind(nlist
, dlist
, &ndih
);
1465 fprintf(stderr
, "%d dihedrals found\n", ndih
);
1469 /* COMPUTE ALL DIHEDRALS! */
1470 read_ang_dih(ftp2fn(efTRX
, NFILE
, fnm
), FALSE
, TRUE
, FALSE
, bPBC
, 1, &idum
,
1471 &nf
, &time
, isize
, index
, &trans_frac
, &aver_angle
, dih
, oenv
);
1473 dt
= (time
[nf
-1]-time
[0])/(nf
-1); /* might want this for corr or n. transit*/
1478 gmx_fatal(FARGS
, "Need at least 2 frames for correlation");
1482 /* put angles in -M_PI to M_PI ! and correct phase factor for phi and psi
1483 * pass nactdih instead of ndih to low_ana_dih_trans and get_chi_product_traj
1484 * to prevent accessing off end of arrays when maxchi < 5 or 6. */
1485 nactdih
= reset_em_all(nlist
, dlist
, nf
, dih
, maxchi
);
1489 dump_em_all(nlist
, dlist
, nf
, time
, dih
, maxchi
, bPhi
, bPsi
, bChi
, bOmega
, bRAD
, oenv
);
1492 /* Histogramming & J coupling constants & calc of S2 order params */
1493 histogramming(log
, nbin
, rt
, nf
, maxchi
, dih
, nlist
, dlist
, index
,
1494 bPhi
, bPsi
, bOmega
, bChi
,
1495 bNormHisto
, bSSHisto
, ftp2fn(efDAT
, NFILE
, fnm
), bfac_max
, &atoms
,
1496 bDo_jc
, opt2fn("-jc", NFILE
, fnm
), oenv
);
1500 * added multiplicity */
1502 snew(multiplicity
, ndih
);
1503 mk_multiplicity_lookup(multiplicity
, maxchi
, nlist
, dlist
, ndih
);
1505 std::strcpy(grpname
, "All residues, ");
1508 std::strcat(grpname
, "Phi ");
1512 std::strcat(grpname
, "Psi ");
1516 std::strcat(grpname
, "Omega ");
1520 std::strcat(grpname
, "Chi 1-");
1521 sprintf(grpname
+ std::strlen(grpname
), "%i", maxchi
);
1525 low_ana_dih_trans(bDo_ot
, opt2fn("-ot", NFILE
, fnm
),
1526 bDo_oh
, opt2fn("-oh", NFILE
, fnm
), maxchi
,
1527 dih
, nlist
, dlist
, nf
, nactdih
, grpname
, multiplicity
,
1528 time
, FALSE
, core_frac
, oenv
);
1530 /* Order parameters */
1531 order_params(log
, opt2fn("-o", NFILE
, fnm
), maxchi
, nlist
, dlist
,
1532 ftp2fn_null(efPDB
, NFILE
, fnm
), bfac_init
,
1533 &atoms
, x
, ePBC
, box
, bPhi
, bPsi
, bChi
, oenv
);
1535 /* Print ramachandran maps! */
1538 do_rama(nf
, nlist
, dlist
, dih
, bViol
, bRamOmega
, oenv
);
1543 do_pp2shifts(log
, nf
, nlist
, dlist
, dih
);
1546 /* rprint S^2, transitions, and rotamer occupancies to log */
1547 traj_t_ns
= 0.001 * (time
[nf
-1]-time
[0]);
1548 pr_dlist(log
, nlist
, dlist
, traj_t_ns
, edPrintST
, bPhi
, bPsi
, bChi
, bOmega
, maxchi
);
1549 pr_dlist(log
, nlist
, dlist
, traj_t_ns
, edPrintRO
, bPhi
, bPsi
, bChi
, bOmega
, maxchi
);
1551 /* transitions to xvg */
1554 print_transitions(opt2fn("-rt", NFILE
, fnm
), maxchi
, nlist
, dlist
, traj_t_ns
, oenv
);
1557 /* chi_product trajectories (ie one "rotamer number" for each residue) */
1558 if (bChiProduct
&& bChi
)
1560 snew(chi_lookup
, nlist
);
1561 for (i
= 0; i
< nlist
; i
++)
1563 snew(chi_lookup
[i
], maxchi
);
1565 mk_chi_lookup(chi_lookup
, maxchi
, nlist
, dlist
);
1567 get_chi_product_traj(dih
, nf
, nactdih
,
1568 maxchi
, dlist
, time
, chi_lookup
, multiplicity
,
1569 FALSE
, bNormHisto
, core_frac
, bAll
,
1570 opt2fn("-cp", NFILE
, fnm
), oenv
);
1572 for (i
= 0; i
< nlist
; i
++)
1574 sfree(chi_lookup
[i
]);
1578 /* Correlation comes last because it messes up the angles */
1581 do_dihcorr(opt2fn("-corr", NFILE
, fnm
), nf
, ndih
, dih
, dt
, nlist
, dlist
, time
,
1582 maxchi
, bPhi
, bPsi
, bChi
, bOmega
, oenv
);
1586 do_view(oenv
, opt2fn("-o", NFILE
, fnm
), "-nxy");
1587 do_view(oenv
, opt2fn("-jc", NFILE
, fnm
), "-nxy");
1590 do_view(oenv
, opt2fn("-corr", NFILE
, fnm
), "-nxy");
1593 gmx_residuetype_destroy(rt
);