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.
53 #include "gromacs/fileio/gmxfio.h"
54 #include "gromacs/fileio/strdb.h"
55 #include "gromacs/topology/invblock.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
61 const char gmx_residuetype_undefined[] = "Other";
63 struct gmx_residuetype
72 static gmx_bool gmx_ask_yesno(gmx_bool bASK)
80 c = toupper(fgetc(stdin));
82 while ((c != 'Y') && (c != 'N'));
92 void write_index(const char *outf, t_blocka *b, char **gnames, gmx_bool bDuplicate, int natoms)
97 out = gmx_fio_fopen(outf, "w");
98 /* fprintf(out,"%5d %5d\n",b->nr,b->nra); */
99 for (i = 0; (i < b->nr); i++)
101 fprintf(out, "[ %s ]\n", gnames[i]);
102 for (k = 0, j = b->index[i]; j < b->index[i+1]; j++, k++)
104 fprintf(out, "%4d ", b->a[j]+1);
113 /* Duplicate copy, useful for computational electrophysiology double-layer setups */
116 fprintf(stderr, "Duplicating the whole system with an atom offset of %d atoms.\n", natoms);
117 for (i = 0; (i < b->nr); i++)
119 fprintf(out, "[ %s_copy ]\n", gnames[i]);
120 for (k = 0, j = b->index[i]; j < b->index[i+1]; j++, k++)
122 fprintf(out, "%4d ", b->a[j]+1 + natoms );
135 void add_grp(t_blocka *b, char ***gnames, int nra, atom_id a[], const char *name)
139 srenew(b->index, b->nr+2);
140 srenew(*gnames, b->nr+1);
141 (*gnames)[b->nr] = strdup(name);
143 srenew(b->a, b->nra+nra);
144 for (i = 0; (i < nra); i++)
146 b->a[b->nra++] = a[i];
149 b->index[b->nr] = b->nra;
152 /* compare index in `a' with group in `b' at `index',
153 when `index'<0 it is relative to end of `b' */
154 static gmx_bool grp_cmp(t_blocka *b, int nra, atom_id a[], int index)
160 index = b->nr-1+index;
164 gmx_fatal(FARGS, "no such index group %d in t_blocka (nr=%d)", index, b->nr);
167 if (nra != b->index[index+1] - b->index[index])
171 for (i = 0; i < nra; i++)
173 if (a[i] != b->a[b->index[index]+i])
182 p_status(const char **restype, int nres, const char **typenames, int ntypes)
189 snew(counter, ntypes);
190 for (i = 0; i < ntypes; i++)
194 for (i = 0; i < nres; i++)
197 for (j = 0; j < ntypes; j++)
199 if (!gmx_strcasecmp(restype[i], typenames[j]))
206 for (i = 0; (i < ntypes); i++)
210 printf("There are: %5d %10s residues\n", counter[i], typenames[i]);
219 mk_aid(t_atoms *atoms, const char ** restype, const char * typestring, int *nra, gmx_bool bMatch)
220 /* Make an array of atom_ids for all atoms with residuetypes matching typestring, or the opposite if bMatch is false */
228 for (i = 0; (i < atoms->nr); i++)
230 res = !gmx_strcasecmp(restype[atoms->atom[i].resind], typestring);
250 static void analyse_other(const char ** restype, t_atoms *atoms,
251 t_blocka *gb, char ***gn, gmx_bool bASK, gmx_bool bVerb)
253 restp_t *restp = NULL;
256 atom_id *other_ndx, *aid, *aaid;
257 int i, j, k, l, resind, naid, naaid, natp, nrestp = 0;
259 for (i = 0; (i < atoms->nres); i++)
261 if (gmx_strcasecmp(restype[i], "Protein") && gmx_strcasecmp(restype[i], "DNA") && gmx_strcasecmp(restype[i], "RNA") && gmx_strcasecmp(restype[i], "Water"))
271 printf("Analysing residues not classified as Protein/DNA/RNA/Water and splitting into groups...\n");
273 snew(other_ndx, atoms->nr);
274 for (k = 0; (k < atoms->nr); k++)
276 resind = atoms->atom[k].resind;
277 rname = *atoms->resinfo[resind].name;
278 if (gmx_strcasecmp(restype[resind], "Protein") && gmx_strcasecmp(restype[resind], "DNA") &&
279 gmx_strcasecmp(restype[resind], "RNA") && gmx_strcasecmp(restype[resind], "Water"))
282 for (l = 0; (l < nrestp); l++)
284 if (strcmp(restp[l].rname, rname) == 0)
291 srenew(restp, nrestp+1);
292 restp[nrestp].rname = strdup(rname);
293 restp[nrestp].bNeg = FALSE;
294 restp[nrestp].gname = strdup(rname);
300 for (i = 0; (i < nrestp); i++)
302 snew(aid, atoms->nr);
304 for (j = 0; (j < atoms->nr); j++)
306 rname = *atoms->resinfo[atoms->atom[j].resind].name;
307 if ((strcmp(restp[i].rname, rname) == 0 && !restp[i].bNeg) ||
308 (strcmp(restp[i].rname, rname) != 0 && restp[i].bNeg))
313 add_grp(gb, gn, naid, aid, restp[i].gname);
316 printf("split %s into atoms (y/n) ? ", restp[i].gname);
318 if (gmx_ask_yesno(bASK))
321 for (k = 0; (k < naid); k++)
323 aname = *atoms->atomname[aid[k]];
324 for (l = 0; (l < natp); l++)
326 if (strcmp(aname, attp[l]) == 0)
333 srenew(attp, ++natp);
334 attp[natp-1] = aname;
339 for (l = 0; (l < natp); l++)
343 for (k = 0; (k < naid); k++)
345 aname = *atoms->atomname[aid[k]];
346 if (strcmp(aname, attp[l]) == 0)
348 aaid[naaid++] = aid[k];
351 add_grp(gb, gn, naaid, aaid, attp[l]);
366 * Cata necessary to construct a single (protein) index group in
369 typedef struct gmx_help_make_index_group
371 /** The set of atom names that will be used to form this index group */
372 const char **defining_atomnames;
373 /** Size of the defining_atomnames array */
374 const int num_defining_atomnames;
375 /** Name of this index group */
376 const char *group_name;
377 /** Whether the above atom names name the atoms in the group, or
378 those not in the group */
379 gmx_bool bTakeComplement;
380 /** The index in wholename gives the first item in the arrays of
381 atomnames that should be tested with 'gmx_strncasecmp' in stead of
382 gmx_strcasecmp, or -1 if all items should be tested with strcasecmp
383 This is comparable to using a '*' wildcard at the end of specific
384 atom names, but that is more involved to implement...
387 /** Only create this index group if it differs from the one specified in compareto,
388 where -1 means to always create this group. */
390 } t_gmx_help_make_index_group;
392 static void analyse_prot(const char ** restype, t_atoms *atoms,
393 t_blocka *gb, char ***gn, gmx_bool bASK, gmx_bool bVerb)
395 /* lists of atomnames to be used in constructing index groups: */
396 static const char *pnoh[] = { "H", "HN" };
397 static const char *pnodum[] = {
398 "MN1", "MN2", "MCB1", "MCB2", "MCG1", "MCG2",
399 "MCD1", "MCD2", "MCE1", "MCE2", "MNZ1", "MNZ2"
401 static const char *calpha[] = { "CA" };
402 static const char *bb[] = { "N", "CA", "C" };
403 static const char *mc[] = { "N", "CA", "C", "O", "O1", "O2", "OC1", "OC2", "OT", "OXT" };
404 static const char *mcb[] = { "N", "CA", "CB", "C", "O", "O1", "O2", "OC1", "OC2", "OT", "OXT" };
405 static const char *mch[] = {
406 "N", "CA", "C", "O", "O1", "O2", "OC1", "OC2", "OT", "OXT",
407 "H1", "H2", "H3", "H", "HN"
410 static const t_gmx_help_make_index_group constructing_data[] =
411 {{ NULL, 0, "Protein", TRUE, -1, -1},
412 { pnoh, asize(pnoh), "Protein-H", TRUE, 0, -1},
413 { calpha, asize(calpha), "C-alpha", FALSE, -1, -1},
414 { bb, asize(bb), "Backbone", FALSE, -1, -1},
415 { mc, asize(mc), "MainChain", FALSE, -1, -1},
416 { mcb, asize(mcb), "MainChain+Cb", FALSE, -1, -1},
417 { mch, asize(mch), "MainChain+H", FALSE, -1, -1},
418 { mch, asize(mch), "SideChain", TRUE, -1, -1},
419 { mch, asize(mch), "SideChain-H", TRUE, 11, -1},
420 { pnodum, asize(pnodum), "Prot-Masses", TRUE, -1, 0}, };
421 const int num_index_groups = asize(constructing_data);
425 int nra, nnpres, npres;
427 char ndx_name[STRLEN], *atnm;
432 printf("Analysing Protein...\n");
434 snew(aid, atoms->nr);
436 /* calculate the number of protein residues */
438 for (i = 0; (i < atoms->nres); i++)
440 if (0 == gmx_strcasecmp(restype[i], "Protein"))
445 /* find matching or complement atoms */
446 for (i = 0; (i < (int)num_index_groups); i++)
449 for (n = 0; (n < atoms->nr); n++)
451 if (0 == gmx_strcasecmp(restype[atoms->atom[n].resind], "Protein"))
454 for (j = 0; (j < constructing_data[i].num_defining_atomnames); j++)
456 /* skip digits at beginning of atomname, e.g. 1H */
457 atnm = *atoms->atomname[n];
458 while (isdigit(atnm[0]))
462 if ( (constructing_data[i].wholename == -1) || (j < constructing_data[i].wholename) )
464 if (0 == gmx_strcasecmp(constructing_data[i].defining_atomnames[j], atnm))
471 if (0 == gmx_strncasecmp(constructing_data[i].defining_atomnames[j], atnm, strlen(constructing_data[i].defining_atomnames[j])))
477 if (constructing_data[i].bTakeComplement != match)
483 /* if we want to add this group always or it differs from previous
485 if (-1 == constructing_data[i].compareto || !grp_cmp(gb, nra, aid, constructing_data[i].compareto-i) )
487 add_grp(gb, gn, nra, aid, constructing_data[i].group_name);
493 for (i = 0; (i < (int)num_index_groups); i++)
495 printf("Split %12s into %5d residues (y/n) ? ", constructing_data[i].group_name, npres);
496 if (gmx_ask_yesno(bASK))
500 for (n = 0; ((atoms->atom[n].resind < npres) && (n < atoms->nr)); )
502 resind = atoms->atom[n].resind;
503 for (; ((atoms->atom[n].resind == resind) && (n < atoms->nr)); n++)
506 for (j = 0; (j < constructing_data[i].num_defining_atomnames); j++)
508 if (0 == gmx_strcasecmp(constructing_data[i].defining_atomnames[j], *atoms->atomname[n]))
513 if (constructing_data[i].bTakeComplement != match)
518 /* copy the residuename to the tail of the groupname */
522 ri = &atoms->resinfo[resind];
523 sprintf(ndx_name, "%s_%s%d%c",
524 constructing_data[i].group_name, *ri->name, ri->nr, ri->ic == ' ' ? '\0' : ri->ic);
525 add_grp(gb, gn, nra, aid, ndx_name);
531 printf("Make group with sidechain and C=O swapped (y/n) ? ");
532 if (gmx_ask_yesno(bASK))
534 /* Make swap sidechain C=O index */
537 for (n = 0; ((atoms->atom[n].resind < npres) && (n < atoms->nr)); )
539 resind = atoms->atom[n].resind;
541 for (; ((atoms->atom[n].resind == resind) && (n < atoms->nr)); n++)
543 if (strcmp("CA", *atoms->atomname[n]) == 0)
549 else if (strcmp("C", *atoms->atomname[n]) == 0)
553 gmx_incons("Atom naming problem");
557 else if (strcmp("O", *atoms->atomname[n]) == 0)
561 gmx_incons("Atom naming problem");
565 else if (strcmp("O1", *atoms->atomname[n]) == 0)
569 gmx_incons("Atom naming problem");
579 /* copy the residuename to the tail of the groupname */
582 add_grp(gb, gn, nra, aid, "SwapSC-CO");
593 /* Return 0 if the name was found, otherwise -1.
594 * p_restype is set to a pointer to the type name, or 'Other' if we did not find it.
597 gmx_residuetype_get_type(gmx_residuetype_t rt, const char * resname, const char ** p_restype)
602 for (i = 0; i < rt->n && rc; i++)
604 rc = gmx_strcasecmp(rt->resname[i], resname);
607 *p_restype = (rc == 0) ? rt->restype[i-1] : gmx_residuetype_undefined;
613 gmx_residuetype_add(gmx_residuetype_t rt, const char *newresname, const char *newrestype)
617 const char * p_oldtype;
619 found = !gmx_residuetype_get_type(rt, newresname, &p_oldtype);
621 if (found && gmx_strcasecmp(p_oldtype, newrestype))
623 fprintf(stderr, "Warning: Residue '%s' already present with type '%s' in database, ignoring new type '%s'.",
624 newresname, p_oldtype, newrestype);
629 srenew(rt->resname, rt->n+1);
630 srenew(rt->restype, rt->n+1);
631 rt->resname[rt->n] = strdup(newresname);
632 rt->restype[rt->n] = strdup(newrestype);
641 gmx_residuetype_init(gmx_residuetype_t *prt)
645 char resname[STRLEN], restype[STRLEN], dum[STRLEN];
648 struct gmx_residuetype *rt;
657 db = libopen("residuetypes.dat");
659 while (get_a_line(db, line, STRLEN))
665 if (sscanf(line, "%s %s %s", resname, restype, dum) != 2)
667 gmx_fatal(FARGS, "Incorrect number of columns (2 expected) for line in residuetypes.dat");
669 gmx_residuetype_add(rt, resname, restype);
681 gmx_residuetype_destroy(gmx_residuetype_t rt)
685 for (i = 0; i < rt->n; i++)
687 sfree(rt->resname[i]);
688 sfree(rt->restype[i]);
698 gmx_residuetype_get_alltypes(gmx_residuetype_t rt,
699 const char *** p_typenames,
704 const char ** my_typename;
710 for (i = 0; i < rt->n; i++)
714 for (j = 0; j < n && !found; j++)
716 found = !gmx_strcasecmp(p, my_typename[j]);
721 srenew(my_typename, n+1);
727 *p_typenames = my_typename;
735 gmx_residuetype_is_protein(gmx_residuetype_t rt, const char *resnm)
740 if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
741 gmx_strcasecmp(p_type, "Protein") == 0)
753 gmx_residuetype_is_dna(gmx_residuetype_t rt, const char *resnm)
758 if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
759 gmx_strcasecmp(p_type, "DNA") == 0)
771 gmx_residuetype_is_rna(gmx_residuetype_t rt, const char *resnm)
776 if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
777 gmx_strcasecmp(p_type, "RNA") == 0)
788 /* Return the size of the arrays */
790 gmx_residuetype_get_size(gmx_residuetype_t rt)
795 /* Search for a residuetype with name resnm within the
796 * gmx_residuetype database. Return the index if found,
800 gmx_residuetype_get_index(gmx_residuetype_t rt, const char *resnm)
805 for (i = 0; i < rt->n && rc; i++)
807 rc = gmx_strcasecmp(rt->resname[i], resnm);
810 return (0 == rc) ? i-1 : -1;
813 /* Return the name of the residuetype with the given index, or
814 * NULL if not found. */
816 gmx_residuetype_get_name(gmx_residuetype_t rt, int index)
818 if (index >= 0 && index < rt->n)
820 return rt->resname[index];
830 void analyse(t_atoms *atoms, t_blocka *gb, char ***gn, gmx_bool bASK, gmx_bool bVerb)
832 gmx_residuetype_t rt = NULL;
835 const char ** restype;
841 const char ** p_typename;
848 printf("Analysing residue names:\n");
850 /* Create system group, every single atom */
851 snew(aid, atoms->nr);
852 for (i = 0; i < atoms->nr; i++)
856 add_grp(gb, gn, atoms->nr, aid, "System");
859 /* For every residue, get a pointer to the residue type name */
860 gmx_residuetype_init(&rt);
863 snew(restype, atoms->nres);
866 for (i = 0; i < atoms->nres; i++)
868 resnm = *atoms->resinfo[i].name;
869 gmx_residuetype_get_type(rt, resnm, &(restype[i]));
871 /* Note that this does not lead to a N*N loop, but N*K, where
872 * K is the number of residue _types_, which is small and independent of N.
875 for (k = 0; k < ntypes && !found; k++)
877 found = !strcmp(restype[i], p_typename[k]);
881 srenew(p_typename, ntypes+1);
882 p_typename[ntypes++] = strdup(restype[i]);
888 p_status(restype, atoms->nres, p_typename, ntypes);
891 for (k = 0; k < ntypes; k++)
893 aid = mk_aid(atoms, restype, p_typename[k], &nra, TRUE);
895 /* Check for special types to do fancy stuff with */
897 if (!gmx_strcasecmp(p_typename[k], "Protein") && nra > 0)
901 analyse_prot(restype, atoms, gb, gn, bASK, bVerb);
903 /* Create a Non-Protein group */
904 aid = mk_aid(atoms, restype, "Protein", &nra, FALSE);
905 if ((nra > 0) && (nra < atoms->nr))
907 add_grp(gb, gn, nra, aid, "non-Protein");
911 else if (!gmx_strcasecmp(p_typename[k], "Water") && nra > 0)
913 add_grp(gb, gn, nra, aid, p_typename[k]);
914 /* Add this group as 'SOL' too, for backward compatibility with older gromacs versions */
915 add_grp(gb, gn, nra, aid, "SOL");
919 /* Solvent, create a negated group too */
920 aid = mk_aid(atoms, restype, "Water", &nra, FALSE);
921 if ((nra > 0) && (nra < atoms->nr))
923 add_grp(gb, gn, nra, aid, "non-Water");
930 add_grp(gb, gn, nra, aid, p_typename[k]);
932 analyse_other(restype, atoms, gb, gn, bASK, bVerb);
938 gmx_residuetype_destroy(rt);
940 /* Create a merged water_and_ions group */
946 for (i = 0; i < gb->nr; i++)
948 if (!gmx_strcasecmp((*gn)[i], "Water"))
951 nwater = gb->index[i+1]-gb->index[i];
953 else if (!gmx_strcasecmp((*gn)[i], "Ion"))
956 nion = gb->index[i+1]-gb->index[i];
960 if (nwater > 0 && nion > 0)
962 srenew(gb->index, gb->nr+2);
963 srenew(*gn, gb->nr+1);
964 (*gn)[gb->nr] = strdup("Water_and_ions");
965 srenew(gb->a, gb->nra+nwater+nion);
968 for (i = gb->index[iwater]; i < gb->index[iwater+1]; i++)
970 gb->a[gb->nra++] = gb->a[i];
975 for (i = gb->index[iion]; i < gb->index[iion+1]; i++)
977 gb->a[gb->nra++] = gb->a[i];
981 gb->index[gb->nr] = gb->nra;
986 void check_index(char *gname, int n, atom_id index[], char *traj, int natoms)
990 for (i = 0; i < n; i++)
992 if (index[i] >= natoms)
994 gmx_fatal(FARGS, "%s atom number (index[%d]=%d) is larger than the number of atoms in %s (%d)",
995 gname ? gname : "Index", i+1, index[i]+1,
996 traj ? traj : "the trajectory", natoms);
998 else if (index[i] < 0)
1000 gmx_fatal(FARGS, "%s atom number (index[%d]=%d) is less than zero",
1001 gname ? gname : "Index", i+1, index[i]+1);
1006 t_blocka *init_index(const char *gfile, char ***grpname)
1011 int i, j, ng, nread;
1012 char line[STRLEN], *pt, str[STRLEN];
1014 in = gmx_fio_fopen(gfile, "r");
1016 get_a_line(in, line, STRLEN);
1028 if (get_header(line, str))
1031 srenew(b->index, b->nr+1);
1032 srenew(*grpname, b->nr);
1037 b->index[b->nr] = b->index[b->nr-1];
1038 (*grpname)[b->nr-1] = strdup(str);
1044 gmx_fatal(FARGS, "The first header of your indexfile is invalid");
1047 while (sscanf(pt, "%s", str) == 1)
1049 i = b->index[b->nr];
1050 if (i >= maxentries)
1053 srenew(b->a, maxentries);
1055 b->a[i] = strtol(str, NULL, 10)-1;
1058 pt = strstr(pt, str)+strlen(str);
1062 while (get_a_line(in, line, STRLEN));
1067 sscanf(line, "%d%d", &b->nr, &b->nra);
1068 snew(b->index, b->nr+1);
1069 snew(*grpname, b->nr);
1072 for (i = 0; (i < b->nr); i++)
1074 nread = fscanf(in, "%s%d", str, &ng);
1075 (*grpname)[i] = strdup(str);
1076 b->index[i+1] = b->index[i]+ng;
1077 if (b->index[i+1] > b->nra)
1079 gmx_fatal(FARGS, "Something wrong in your indexfile at group %s", str);
1081 for (j = 0; (j < ng); j++)
1083 nread = fscanf(in, "%d", &a);
1084 b->a[b->index[i]+j] = a;
1090 for (i = 0; (i < b->nr); i++)
1092 for (j = b->index[i]; (j < b->index[i+1]); j++)
1096 fprintf(stderr, "\nWARNING: negative index %d in group %s\n\n",
1097 b->a[j], (*grpname)[i]);
1105 static void minstring(char *str)
1109 for (i = 0; (i < (int)strlen(str)); i++)
1118 int find_group(char s[], int ngrps, char **grpname)
1121 char string[STRLEN];
1127 /* first look for whole name match */
1130 for (i = 0; i < ngrps; i++)
1132 if (gmx_strcasecmp_min(s, grpname[i]) == 0)
1142 /* second look for first string match */
1145 for (i = 0; i < ngrps; i++)
1147 if (gmx_strncasecmp_min(s, grpname[i], n) == 0)
1157 /* last look for arbitrary substring match */
1162 for (i = 0; i < ngrps; i++)
1164 strcpy(string, grpname[i]);
1167 if (strstr(string, s) != NULL)
1179 printf("Error: Multiple groups '%s' selected\n", s);
1185 static int qgroup(int *a, int ngrps, char **grpname)
1194 fprintf(stderr, "Select a group: ");
1197 if (scanf("%s", s) != 1)
1199 gmx_fatal(FARGS, "Cannot read from input");
1201 trim(s); /* remove spaces */
1203 while (strlen(s) == 0);
1204 aa = strtol(s, &end, 10);
1205 if (aa == 0 && end[0] != '\0') /* string entered */
1207 aa = find_group(s, ngrps, grpname);
1209 bInRange = (aa >= 0 && aa < ngrps);
1212 printf("Error: No such group '%s'\n", s);
1216 printf("Selected %d: '%s'\n", aa, grpname[aa]);
1221 static void rd_groups(t_blocka *grps, char **grpname, char *gnames[],
1222 int ngrps, int isize[], atom_id *index[], int grpnr[])
1228 gmx_fatal(FARGS, "Error: no groups in indexfile");
1230 for (i = 0; (i < grps->nr); i++)
1232 fprintf(stderr, "Group %5d (%15s) has %5d elements\n", i, grpname[i],
1233 grps->index[i+1]-grps->index[i]);
1235 for (i = 0; (i < ngrps); i++)
1241 gnr1 = qgroup(&grpnr[i], grps->nr, grpname);
1242 if ((gnr1 < 0) || (gnr1 >= grps->nr))
1244 fprintf(stderr, "Select between %d and %d.\n", 0, grps->nr-1);
1247 while ((gnr1 < 0) || (gnr1 >= grps->nr));
1251 fprintf(stderr, "There is one group in the index\n");
1254 gnames[i] = strdup(grpname[gnr1]);
1255 isize[i] = grps->index[gnr1+1]-grps->index[gnr1];
1256 snew(index[i], isize[i]);
1257 for (j = 0; (j < isize[i]); j++)
1259 index[i][j] = grps->a[grps->index[gnr1]+j];
1264 void rd_index(const char *statfile, int ngrps, int isize[],
1265 atom_id *index[], char *grpnames[])
1274 gmx_fatal(FARGS, "No index file specified");
1276 grps = init_index(statfile, &gnames);
1277 rd_groups(grps, gnames, grpnames, ngrps, isize, index, grpnr);
1280 void rd_index_nrs(char *statfile, int ngrps, int isize[],
1281 atom_id *index[], char *grpnames[], int grpnr[])
1288 gmx_fatal(FARGS, "No index file specified");
1290 grps = init_index(statfile, &gnames);
1292 rd_groups(grps, gnames, grpnames, ngrps, isize, index, grpnr);
1295 void get_index(t_atoms *atoms, const char *fnm, int ngrps,
1296 int isize[], atom_id *index[], char *grpnames[])
1299 t_blocka *grps = NULL;
1306 grps = init_index(fnm, gnames);
1311 snew(grps->index, 1);
1312 analyse(atoms, grps, gnames, FALSE, FALSE);
1316 gmx_incons("You need to supply a valid atoms structure or a valid index file name");
1319 rd_groups(grps, *gnames, grpnames, ngrps, isize, index, grpnr);
1322 t_cluster_ndx *cluster_index(FILE *fplog, const char *ndx)
1328 c->clust = init_index(ndx, &c->grpname);
1330 for (i = 0; (i < c->clust->nra); i++)
1332 c->maxframe = max(c->maxframe, c->clust->a[i]);
1334 fprintf(fplog ? fplog : stdout,
1335 "There are %d clusters containing %d structures, highest framenr is %d\n",
1336 c->clust->nr, c->clust->nra, c->maxframe);
1339 pr_blocka(debug, 0, "clust", c->clust, TRUE);
1340 for (i = 0; (i < c->clust->nra); i++)
1342 if ((c->clust->a[i] < 0) || (c->clust->a[i] > c->maxframe))
1344 gmx_fatal(FARGS, "Range check error for c->clust->a[%d] = %d\n"
1345 "should be within 0 and %d", i, c->clust->a[i], c->maxframe+1);
1349 c->inv_clust = make_invblocka(c->clust, c->maxframe);