Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / gmxana / dlist.cpp
blob70ded720444045a5e5fc918e5b4db3795c23ef41
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,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 <cstdlib>
40 #include <cstring>
42 #include "gromacs/gmxana/gstat.h"
43 #include "gromacs/topology/residuetypes.h"
44 #include "gromacs/topology/topology.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/smalloc.h"
48 t_dlist *mk_dlist(FILE *log,
49 const t_atoms *atoms, int *nlist,
50 gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bHChi,
51 int maxchi, int r0, ResidueType *rt)
53 int i, j, ii;
54 t_dihatms atm, prev;
55 int nl = 0, nc[edMax];
56 char *thisres;
57 t_dlist *dl;
59 snew(dl, atoms->nres+1);
60 prev.C = prev.Cn[1] = -1; /* Keep the compiler quiet */
61 for (i = 0; (i < edMax); i++)
63 nc[i] = 0;
65 i = 0;
66 while (i < atoms->nr)
68 int 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++)
74 atm.Cn[j] = -1;
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 ((std::strcmp(*(atoms->atomname[i]), "H") == 0) ||
82 (std::strcmp(*(atoms->atomname[i]), "H1") == 0) ||
83 (std::strcmp(*(atoms->atomname[i]), "HN") == 0) )
85 atm.H = i;
87 else if (std::strcmp(*(atoms->atomname[i]), "N") == 0)
89 atm.N = i;
91 else if (std::strcmp(*(atoms->atomname[i]), "C") == 0)
93 atm.C = i;
95 else if ((std::strcmp(*(atoms->atomname[i]), "O") == 0) ||
96 (std::strcmp(*(atoms->atomname[i]), "O1") == 0) ||
97 (std::strcmp(*(atoms->atomname[i]), "OC1") == 0) ||
98 (std::strcmp(*(atoms->atomname[i]), "OT1") == 0))
100 atm.O = i;
102 else if (std::strcmp(*(atoms->atomname[i]), "CA") == 0)
104 atm.Cn[1] = i;
106 else if (std::strcmp(*(atoms->atomname[i]), "CB") == 0)
108 atm.Cn[2] = i;
110 else if ((std::strcmp(*(atoms->atomname[i]), "CG") == 0) ||
111 (std::strcmp(*(atoms->atomname[i]), "CG1") == 0) ||
112 (std::strcmp(*(atoms->atomname[i]), "OG") == 0) ||
113 (std::strcmp(*(atoms->atomname[i]), "OG1") == 0) ||
114 (std::strcmp(*(atoms->atomname[i]), "SG") == 0))
116 atm.Cn[3] = i;
118 else if ((std::strcmp(*(atoms->atomname[i]), "CD") == 0) ||
119 (std::strcmp(*(atoms->atomname[i]), "CD1") == 0) ||
120 (std::strcmp(*(atoms->atomname[i]), "SD") == 0) ||
121 (std::strcmp(*(atoms->atomname[i]), "OD1") == 0) ||
122 (std::strcmp(*(atoms->atomname[i]), "ND1") == 0))
124 atm.Cn[4] = i;
126 /* by grs - split the Cn[4] into 2 bits to check allowing dih to H */
127 else if (bHChi && ((std::strcmp(*(atoms->atomname[i]), "HG") == 0) ||
128 (std::strcmp(*(atoms->atomname[i]), "HG1") == 0)) )
130 atm.Cn[4] = i;
132 else if ((std::strcmp(*(atoms->atomname[i]), "CE") == 0) ||
133 (std::strcmp(*(atoms->atomname[i]), "CE1") == 0) ||
134 (std::strcmp(*(atoms->atomname[i]), "OE1") == 0) ||
135 (std::strcmp(*(atoms->atomname[i]), "NE") == 0))
137 atm.Cn[5] = i;
139 else if ((std::strcmp(*(atoms->atomname[i]), "CZ") == 0) ||
140 (std::strcmp(*(atoms->atomname[i]), "NZ") == 0))
142 atm.Cn[6] = i;
144 /* HChi flag here too */
145 else if (bHChi && (std::strcmp(*(atoms->atomname[i]), "NH1") == 0))
147 atm.Cn[7] = i;
149 i++;
152 thisres = *(atoms->resinfo[ires].name);
154 /* added by grs - special case for aromatics, whose chis above 2 are
155 not real and produce rubbish output - so set back to -1 */
156 if (std::strcmp(thisres, "PHE") == 0 ||
157 std::strcmp(thisres, "TYR") == 0 ||
158 std::strcmp(thisres, "PTR") == 0 ||
159 std::strcmp(thisres, "TRP") == 0 ||
160 std::strcmp(thisres, "HIS") == 0 ||
161 std::strcmp(thisres, "HISA") == 0 ||
162 std::strcmp(thisres, "HISB") == 0)
164 for (ii = 5; ii <= 7; ii++)
166 atm.Cn[ii] = -1;
169 /* end fixing aromatics */
171 /* Special case for Pro, has no H */
172 if (std::strcmp(thisres, "PRO") == 0)
174 atm.H = atm.Cn[4];
176 /* Carbon from previous residue */
177 if (prev.C != -1)
179 atm.minC = prev.C;
181 /* Alpha-carbon from previous residue */
182 if (prev.Cn[1] != -1)
184 atm.minCalpha = prev.Cn[1];
186 prev = atm;
188 /* Check how many dihedrals we have */
189 if ((atm.N != -1) && (atm.Cn[1] != -1) && (atm.C != -1) &&
190 (atm.O != -1) && ((atm.H != -1) || (atm.minC != -1)))
192 dl[nl].resnr = ires+1;
193 dl[nl].atm = atm;
194 dl[nl].atm.Cn[0] = atm.N;
195 if ((atm.Cn[3] != -1) && (atm.Cn[2] != -1) && (atm.Cn[1] != -1))
197 nc[0]++;
198 if (atm.Cn[4] != -1)
200 nc[1]++;
201 if (atm.Cn[5] != -1)
203 nc[2]++;
204 if (atm.Cn[6] != -1)
206 nc[3]++;
207 if (atm.Cn[7] != -1)
209 nc[4]++;
210 if (atm.Cn[8] != -1)
212 nc[5]++;
219 if ((atm.minC != -1) && (atm.minCalpha != -1))
221 nc[6]++;
223 dl[nl].index = rt->indexFromResidueName(thisres);
225 /* Prevent seg fault from unknown residues. If one adds a custom residue to
226 * residuetypes.dat but somehow loses it, changes it, or does analysis on
227 * another machine, the residue type will be unknown. */
228 if (dl[nl].index == -1)
230 gmx_fatal(FARGS, "Unknown residue %s when searching for residue type.\n"
231 "Maybe you need to add a custom residue in residuetypes.dat.",
232 thisres);
235 sprintf(dl[nl].name, "%s%d", thisres, ires+r0);
236 nl++;
238 else if (debug)
240 fprintf(debug, "Could not find N atom but could find other atoms"
241 " in residue %s%d\n", thisres, ires+r0);
244 fprintf(stderr, "\n");
245 fprintf(log, "\n");
246 fprintf(log, "There are %d residues with dihedrals\n", nl);
247 j = 0;
248 if (bPhi)
250 j += nl;
252 if (bPsi)
254 j += nl;
256 if (bChi)
258 for (i = 0; (i < maxchi); i++)
260 j += nc[i];
263 fprintf(log, "There are %d dihedrals\n", j);
264 fprintf(log, "Dihedral: ");
265 if (bPhi)
267 fprintf(log, " Phi ");
269 if (bPsi)
271 fprintf(log, " Psi ");
273 if (bChi)
275 for (i = 0; (i < maxchi); i++)
277 fprintf(log, "Chi%d ", i+1);
280 fprintf(log, "\nNumber: ");
281 if (bPhi)
283 fprintf(log, "%4d ", nl);
285 if (bPsi)
287 fprintf(log, "%4d ", nl);
289 if (bChi)
291 for (i = 0; (i < maxchi); i++)
293 fprintf(log, "%4d ", nc[i]);
296 fprintf(log, "\n");
298 *nlist = nl;
300 return dl;
303 gmx_bool has_dihedral(int Dih, t_dlist *dl)
305 gmx_bool b = FALSE;
306 int ddd;
308 switch (Dih)
310 case edPhi:
311 b = ((dl->atm.H != -1) && (dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1));
312 break;
313 case edPsi:
314 b = ((dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1) && (dl->atm.O != -1));
315 break;
316 case edOmega:
317 b = ((dl->atm.minCalpha != -1) && (dl->atm.minC != -1) && (dl->atm.N != -1) && (dl->atm.Cn[1] != -1));
318 break;
319 case edChi1:
320 case edChi2:
321 case edChi3:
322 case edChi4:
323 case edChi5:
324 case edChi6:
325 ddd = Dih - edChi1;
326 b = ((dl->atm.Cn[ddd] != -1) && (dl->atm.Cn[ddd+1] != -1) &&
327 (dl->atm.Cn[ddd+2] != -1) && (dl->atm.Cn[ddd+3] != -1));
328 break;
329 default:
330 pr_dlist(stdout, 1, dl, 1, 0, TRUE, TRUE, TRUE, TRUE, MAXCHI);
331 gmx_fatal(FARGS, "Non existent dihedral %d in file %s, line %d",
332 Dih, __FILE__, __LINE__);
334 return b;
337 static void pr_one_ro(FILE *fp, t_dlist *dl, int nDih, real gmx_unused dt)
339 int k;
340 for (k = 0; k < NROT; k++)
342 fprintf(fp, " %6.2f", dl->rot_occ[nDih][k]);
344 fprintf(fp, "\n");
347 static void pr_ntr_s2(FILE *fp, t_dlist *dl, int nDih, real dt)
349 fprintf(fp, " %6.2f %6.2f\n", (dt == 0) ? 0 : dl->ntr[nDih]/dt, dl->S2[nDih]);
352 void pr_dlist(FILE *fp, int nl, t_dlist dl[], real dt, int printtype,
353 gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bOmega, int maxchi)
355 int i, Xi;
357 void (*pr_props)(FILE *, t_dlist *, int, real);
359 /* Analysis of dihedral transitions etc */
361 if (printtype == edPrintST)
363 pr_props = pr_ntr_s2;
364 fprintf(stderr, "Now printing out transitions and OPs...\n");
366 else
368 pr_props = pr_one_ro;
369 fprintf(stderr, "Now printing out rotamer occupancies...\n");
370 fprintf(fp, "\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX\n\n");
373 /* change atom numbers from 0 based to 1 based */
374 for (i = 0; (i < nl); i++)
376 fprintf(fp, "Residue %s\n", dl[i].name);
377 if (printtype == edPrintST)
379 fprintf(fp, " Angle [ AI, AJ, AK, AL] #tr/ns S^2D \n"
380 "--------------------------------------------\n");
382 else
384 fprintf(fp, " Angle [ AI, AJ, AK, AL] rotamers 0 g(-) t g(+)\n"
385 "--------------------------------------------\n");
387 if (bPhi)
389 fprintf(fp, " Phi [%5d,%5d,%5d,%5d]",
390 (dl[i].atm.H == -1) ? 1+dl[i].atm.minC : 1+dl[i].atm.H,
391 1+dl[i].atm.N, 1+dl[i].atm.Cn[1], 1+dl[i].atm.C);
392 pr_props(fp, &dl[i], edPhi, dt);
394 if (bPsi)
396 fprintf(fp, " Psi [%5d,%5d,%5d,%5d]", 1+dl[i].atm.N, 1+dl[i].atm.Cn[1],
397 1+dl[i].atm.C, 1+dl[i].atm.O);
398 pr_props(fp, &dl[i], edPsi, dt);
400 if (bOmega && has_dihedral(edOmega, &(dl[i])))
402 fprintf(fp, " Omega [%5d,%5d,%5d,%5d]", 1+dl[i].atm.minCalpha, 1+dl[i].atm.minC,
403 1+dl[i].atm.N, 1+dl[i].atm.Cn[1]);
404 pr_props(fp, &dl[i], edOmega, dt);
406 for (Xi = 0; Xi < MAXCHI; Xi++)
408 if (bChi && (Xi < maxchi) && (dl[i].atm.Cn[Xi+3] != -1) )
410 fprintf(fp, " Chi%d[%5d,%5d,%5d,%5d]", Xi+1, 1+dl[i].atm.Cn[Xi],
411 1+dl[i].atm.Cn[Xi+1], 1+dl[i].atm.Cn[Xi+2],
412 1+dl[i].atm.Cn[Xi+3]);
413 pr_props(fp, &dl[i], Xi+edChi1, dt); /* Xi+2 was wrong here */
416 fprintf(fp, "\n");
422 int pr_trans(FILE *fp, int nl, t_dlist dl[], real dt, int Xi)
424 /* never called at the moment */
426 int i, nn, nz;
428 nz = 0;
429 fprintf(fp, "\\begin{table}[h]\n");
430 fprintf(fp, "\\caption{Number of dihedral transitions per nanosecond}\n");
431 fprintf(fp, "\\begin{tabular}{|l|l|}\n");
432 fprintf(fp, "\\hline\n");
433 fprintf(fp, "Residue\t&$\\chi_%d$\t\\\\\n", Xi+1);
434 for (i = 0; (i < nl); i++)
436 nn = static_cast<int>(dl[i].ntr[Xi]/dt);
438 if (nn == 0)
440 fprintf(fp, "%s\t&\\HL{%d}\t\\\\\n", dl[i].name, nn);
441 nz++;
443 else if (nn > 0)
445 fprintf(fp, "%s\t&\\%d\t\\\\\n", dl[i].name, nn);
448 fprintf(fp, "\\hline\n");
449 fprintf(fp, "\\end{tabular}\n");
450 fprintf(fp, "\\end{table}\n\n");
452 return nz;