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 * GROningen Mixture of Alchemy and Childrens' Stories
62 const char gmx_residuetype_undefined[]="Other";
64 struct gmx_residuetype
73 static gmx_bool gmx_ask_yesno(gmx_bool bASK)
79 c=toupper(fgetc(stdin));
80 } while ((c != 'Y') && (c != 'N'));
88 t_blocka *new_blocka(void)
98 void write_index(const char *outf, t_blocka *b,char **gnames)
103 out=gmx_fio_fopen(outf,"w");
104 /* fprintf(out,"%5d %5d\n",b->nr,b->nra); */
105 for(i=0; (i<b->nr); i++) {
106 fprintf(out,"[ %s ]\n",gnames[i]);
107 for(k=0,j=b->index[i]; j<b->index[i+1]; j++,k++) {
108 fprintf(out,"%4d ",b->a[j]+1);
117 void add_grp(t_blocka *b,char ***gnames,int nra,atom_id a[],const char *name)
121 srenew(b->index,b->nr+2);
122 srenew(*gnames,b->nr+1);
123 (*gnames)[b->nr]=strdup(name);
125 srenew(b->a,b->nra+nra);
126 for(i=0; (i<nra); i++)
129 b->index[b->nr]=b->nra;
132 /* compare index in `a' with group in `b' at `index',
133 when `index'<0 it is relative to end of `b' */
134 static gmx_bool grp_cmp(t_blocka *b, int nra, atom_id a[], int index)
139 index = b->nr-1+index;
141 gmx_fatal(FARGS,"no such index group %d in t_blocka (nr=%d)",index,b->nr);
143 if ( nra != b->index[index+1] - b->index[index] )
146 if ( a[i] != b->a[b->index[index]+i] )
152 p_status(const char **restype, int nres, const char **typenames, int ntypes)
159 snew(counter,ntypes);
160 for(i=0;i<ntypes;i++)
167 for(j=0;j<ntypes;j++)
169 if(!gmx_strcasecmp(restype[i],typenames[j]))
176 for(i=0; (i<ntypes); i++)
180 printf("There are: %5d %10s residues\n",counter[i],typenames[i]);
189 mk_aid(t_atoms *atoms,const char ** restype,const char * typestring,int *nra,gmx_bool bMatch)
190 /* Make an array of atom_ids for all atoms with residuetypes matching typestring, or the opposite if bMatch is false */
198 for(i=0; (i<atoms->nr); i++)
200 res=!gmx_strcasecmp(restype[atoms->atom[i].resind],typestring);
220 static void analyse_other(const char ** restype,t_atoms *atoms,
221 t_blocka *gb,char ***gn,gmx_bool bASK,gmx_bool bVerb)
226 atom_id *other_ndx,*aid,*aaid;
227 int i,j,k,l,resind,naid,naaid,natp,nrestp=0;
229 for(i=0; (i<atoms->nres); i++)
231 if (gmx_strcasecmp(restype[i],"Protein") && gmx_strcasecmp(restype[i],"DNA") && gmx_strcasecmp(restype[i],"RNA") && gmx_strcasecmp(restype[i],"Water"))
236 if (i < atoms->nres) {
239 printf("Analysing residues not classified as Protein/DNA/RNA/Water and splitting into groups...\n");
240 snew(other_ndx,atoms->nr);
241 for(k=0; (k<atoms->nr); k++) {
242 resind = atoms->atom[k].resind;
243 rname = *atoms->resinfo[resind].name;
244 if (gmx_strcasecmp(restype[resind],"Protein") && gmx_strcasecmp(restype[resind],"DNA") &&
245 gmx_strcasecmp(restype[resind],"RNA") && gmx_strcasecmp(restype[resind],"Water"))
248 for(l=0; (l<nrestp); l++)
249 if (strcmp(restp[l].rname,rname) == 0)
252 srenew(restp,nrestp+1);
253 restp[nrestp].rname = strdup(rname);
254 restp[nrestp].bNeg = FALSE;
255 restp[nrestp].gname = strdup(rname);
261 for(i=0; (i<nrestp); i++) {
264 for(j=0; (j<atoms->nr); j++) {
265 rname = *atoms->resinfo[atoms->atom[j].resind].name;
266 if ((strcmp(restp[i].rname,rname) == 0 && !restp[i].bNeg) ||
267 (strcmp(restp[i].rname,rname) != 0 && restp[i].bNeg)) {
271 add_grp(gb,gn,naid,aid,restp[i].gname);
273 printf("split %s into atoms (y/n) ? ",restp[i].gname);
275 if (gmx_ask_yesno(bASK)) {
277 for(k=0; (k<naid); k++) {
278 aname=*atoms->atomname[aid[k]];
279 for(l=0; (l<natp); l++)
280 if (strcmp(aname,attp[l]) == 0)
288 for(l=0; (l<natp); l++) {
291 for(k=0; (k<naid); k++) {
292 aname=*atoms->atomname[aid[k]];
293 if (strcmp(aname,attp[l])==0)
294 aaid[naaid++]=aid[k];
296 add_grp(gb,gn,naaid,aaid,attp[l]);
310 /*! /brief Instances of this struct contain the data necessary to
311 * construct a single (protein) index group in
313 typedef struct gmx_help_make_index_group
315 /* The set of atom names that will be used to form this index group */
316 const char **defining_atomnames;
317 /* Size of the defining_atomnames array */
318 const int num_defining_atomnames;
319 /* Name of this index group */
320 const char *group_name;
321 /* Whether the above atom names name the atoms in the group, or
322 * those not in the group */
323 gmx_bool bTakeComplement;
324 /* The index in wholename gives the first item in the arrays of
325 * atomnames that should be tested with 'gmx_strncasecmp' in stead of
326 * gmx_strcasecmp, or -1 if all items should be tested with strcasecmp
327 * This is comparable to using a '*' wildcard at the end of specific
328 * atom names, but that is more involved to implement...
331 /* Only create this index group if it differs from the one specified in compareto,
332 where -1 means to always create this group. */
334 } t_gmx_help_make_index_group;
336 static void analyse_prot(const char ** restype,t_atoms *atoms,
337 t_blocka *gb,char ***gn,gmx_bool bASK,gmx_bool bVerb)
339 /* lists of atomnames to be used in constructing index groups: */
340 static const char *pnoh[] = { "H", "HN" };
341 static const char *pnodum[] = { "MN1", "MN2", "MCB1", "MCB2", "MCG1", "MCG2",
342 "MCD1", "MCD2", "MCE1", "MCE2", "MNZ1", "MNZ2" };
343 static const char *calpha[] = { "CA" };
344 static const char *bb[] = { "N","CA","C" };
345 static const char *mc[] = { "N","CA","C","O","O1","O2","OC1","OC2","OT","OXT" };
346 static const char *mcb[] = { "N","CA","CB","C","O","O1","O2","OC1","OC2","OT","OXT" };
347 static const char *mch[] = { "N","CA","C","O","O1","O2","OC1","OC2","OT","OXT",
348 "H1","H2","H3","H","HN" };
350 static const t_gmx_help_make_index_group constructing_data[] =
351 {{ NULL, 0, "Protein", TRUE, -1, -1},
352 { pnoh, asize(pnoh), "Protein-H", TRUE, 0, -1},
353 { calpha, asize(calpha), "C-alpha", FALSE, -1, -1},
354 { bb, asize(bb), "Backbone", FALSE, -1, -1},
355 { mc, asize(mc), "MainChain", FALSE, -1, -1},
356 { mcb, asize(mcb), "MainChain+Cb", FALSE, -1, -1},
357 { mch, asize(mch), "MainChain+H", FALSE, -1, -1},
358 { mch, asize(mch), "SideChain", TRUE, -1, -1},
359 { mch, asize(mch), "SideChain-H", TRUE, 11, -1},
360 { pnodum, asize(pnodum), "Prot-Masses", TRUE, -1, 0},
362 const int num_index_groups = asize(constructing_data);
366 int nra,nnpres,npres;
368 char ndx_name[STRLEN],*atnm;
373 printf("Analysing Protein...\n");
377 /* calculate the number of protein residues */
379 for(i=0; (i<atoms->nres); i++) {
380 if (0 == gmx_strcasecmp(restype[i],"Protein")) {
384 /* find matching or complement atoms */
385 for(i=0; (i<(int)num_index_groups); i++) {
387 for(n=0; (n<atoms->nr); n++) {
388 if (0 == gmx_strcasecmp(restype[atoms->atom[n].resind],"Protein")) {
390 for(j=0; (j<constructing_data[i].num_defining_atomnames); j++) {
391 /* skip digits at beginning of atomname, e.g. 1H */
392 atnm=*atoms->atomname[n];
393 while (isdigit(atnm[0])) {
396 if ( (constructing_data[i].wholename==-1) || (j<constructing_data[i].wholename) ) {
397 if (0 == gmx_strcasecmp(constructing_data[i].defining_atomnames[j],atnm)) {
401 if (0 == gmx_strncasecmp(constructing_data[i].defining_atomnames[j],atnm,strlen(constructing_data[i].defining_atomnames[j]))) {
406 if (constructing_data[i].bTakeComplement != match) {
411 /* if we want to add this group always or it differs from previous
413 if ( -1 == constructing_data[i].compareto || !grp_cmp(gb,nra,aid,constructing_data[i].compareto-i) ) {
414 add_grp(gb,gn,nra,aid,constructing_data[i].group_name);
419 for(i=0; (i<(int)num_index_groups); i++) {
420 printf("Split %12s into %5d residues (y/n) ? ",constructing_data[i].group_name,npres);
421 if (gmx_ask_yesno(bASK)) {
424 for(n=0;((atoms->atom[n].resind < npres) && (n<atoms->nr));) {
425 resind = atoms->atom[n].resind;
426 for(;((atoms->atom[n].resind==resind) && (n<atoms->nr));n++) {
428 for(j=0;(j<constructing_data[i].num_defining_atomnames); j++) {
429 if (0 == gmx_strcasecmp(constructing_data[i].defining_atomnames[j],*atoms->atomname[n])) {
433 if (constructing_data[i].bTakeComplement != match) {
437 /* copy the residuename to the tail of the groupname */
440 ri = &atoms->resinfo[resind];
441 sprintf(ndx_name,"%s_%s%d%c",
442 constructing_data[i].group_name,*ri->name,ri->nr,ri->ic==' ' ? '\0' : ri->ic);
443 add_grp(gb,gn,nra,aid,ndx_name);
449 printf("Make group with sidechain and C=O swapped (y/n) ? ");
450 if (gmx_ask_yesno(bASK)) {
451 /* Make swap sidechain C=O index */
454 for(n=0;((atoms->atom[n].resind < npres) && (n<atoms->nr));) {
455 resind = atoms->atom[n].resind;
457 for(;((atoms->atom[n].resind==resind) && (n<atoms->nr));n++)
458 if (strcmp("CA",*atoms->atomname[n]) == 0) {
462 } else if (strcmp("C",*atoms->atomname[n]) == 0) {
464 gmx_incons("Atom naming problem");
467 } else if (strcmp("O",*atoms->atomname[n]) == 0) {
469 gmx_incons("Atom naming problem");
472 } else if (strcmp("O1",*atoms->atomname[n]) == 0) {
474 gmx_incons("Atom naming problem");
480 /* copy the residuename to the tail of the groupname */
482 add_grp(gb,gn,nra,aid,"SwapSC-CO");
493 /* Return 0 if the name was found, otherwise -1.
494 * p_restype is set to a pointer to the type name, or 'Other' if we did not find it.
497 gmx_residuetype_get_type(gmx_residuetype_t rt,const char * resname, const char ** p_restype)
502 for(i=0;i<rt->n && rc;i++)
504 rc=gmx_strcasecmp(rt->resname[i],resname);
507 *p_restype = (rc==0) ? rt->restype[i-1] : gmx_residuetype_undefined;
513 gmx_residuetype_add(gmx_residuetype_t rt,const char *newresname, const char *newrestype)
517 const char * p_oldtype;
519 found = !gmx_residuetype_get_type(rt,newresname,&p_oldtype);
521 if(found && gmx_strcasecmp(p_oldtype,newrestype))
523 fprintf(stderr,"Warning: Residue '%s' already present with type '%s' in database, ignoring new type '%s'.",
524 newresname,p_oldtype,newrestype);
529 srenew(rt->resname,rt->n+1);
530 srenew(rt->restype,rt->n+1);
531 rt->resname[rt->n]=strdup(newresname);
532 rt->restype[rt->n]=strdup(newrestype);
541 gmx_residuetype_init(gmx_residuetype_t *prt)
545 char resname[STRLEN],restype[STRLEN],dum[STRLEN];
548 struct gmx_residuetype *rt;
557 db=libopen("residuetypes.dat");
559 while(get_a_line(db,line,STRLEN))
565 if(sscanf(line,"%s %s %s",resname,restype,dum)!=2)
567 gmx_fatal(FARGS,"Incorrect number of columns (2 expected) for line in residuetypes.dat");
569 gmx_residuetype_add(rt,resname,restype);
581 gmx_residuetype_destroy(gmx_residuetype_t rt)
587 sfree(rt->resname[i]);
588 sfree(rt->restype[i]);
598 gmx_residuetype_get_alltypes(gmx_residuetype_t rt,
599 const char *** p_typenames,
604 const char ** my_typename;
614 for(j=0;j<n && !found;j++)
616 found=!gmx_strcasecmp(p,my_typename[j]);
621 srenew(my_typename,n+1);
627 *p_typenames=my_typename;
635 gmx_residuetype_is_protein(gmx_residuetype_t rt, const char *resnm)
640 if(gmx_residuetype_get_type(rt,resnm,&p_type)==0 &&
641 gmx_strcasecmp(p_type,"Protein")==0)
653 gmx_residuetype_is_dna(gmx_residuetype_t rt, const char *resnm)
658 if(gmx_residuetype_get_type(rt,resnm,&p_type)==0 &&
659 gmx_strcasecmp(p_type,"DNA")==0)
671 gmx_residuetype_is_rna(gmx_residuetype_t rt, const char *resnm)
676 if(gmx_residuetype_get_type(rt,resnm,&p_type)==0 &&
677 gmx_strcasecmp(p_type,"RNA")==0)
688 /* Return the size of the arrays */
690 gmx_residuetype_get_size(gmx_residuetype_t rt)
695 /* Search for a residuetype with name resnm within the
696 * gmx_residuetype database. Return the index if found,
700 gmx_residuetype_get_index(gmx_residuetype_t rt, const char *resnm)
705 for(i=0;i<rt->n && rc;i++)
707 rc=gmx_strcasecmp(rt->resname[i],resnm);
710 return (0 == rc) ? i-1 : -1;
713 /* Return the name of the residuetype with the given index, or
714 * NULL if not found. */
716 gmx_residuetype_get_name(gmx_residuetype_t rt, int index)
718 if(index >= 0 && index < rt->n) {
719 return rt->resname[index];
727 void analyse(t_atoms *atoms,t_blocka *gb,char ***gn,gmx_bool bASK,gmx_bool bVerb)
729 gmx_residuetype_t rt=NULL;
732 const char ** restype;
738 const char ** p_typename;
745 printf("Analysing residue names:\n");
747 /* Create system group, every single atom */
749 for(i=0;i<atoms->nr;i++)
753 add_grp(gb,gn,atoms->nr,aid,"System");
756 /* For every residue, get a pointer to the residue type name */
757 gmx_residuetype_init(&rt);
760 snew(restype,atoms->nres);
763 for(i=0;i<atoms->nres;i++)
765 resnm = *atoms->resinfo[i].name;
766 gmx_residuetype_get_type(rt,resnm,&(restype[i]));
768 /* Note that this does not lead to a N*N loop, but N*K, where
769 * K is the number of residue _types_, which is small and independent of N.
772 for(k=0;k<ntypes && !found;k++)
774 found = !strcmp(restype[i],p_typename[k]);
778 srenew(p_typename,ntypes+1);
779 p_typename[ntypes++] = strdup(restype[i]);
785 p_status(restype,atoms->nres,p_typename,ntypes);
788 for(k=0;k<ntypes;k++)
790 aid=mk_aid(atoms,restype,p_typename[k],&nra,TRUE);
792 /* Check for special types to do fancy stuff with */
794 if(!gmx_strcasecmp(p_typename[k],"Protein") && nra>0)
798 analyse_prot(restype,atoms,gb,gn,bASK,bVerb);
800 /* Create a Non-Protein group */
801 aid=mk_aid(atoms,restype,"Protein",&nra,FALSE);
802 if ((nra > 0) && (nra < atoms->nr))
804 add_grp(gb,gn,nra,aid,"non-Protein");
808 else if(!gmx_strcasecmp(p_typename[k],"Water") && nra>0)
810 add_grp(gb,gn,nra,aid,p_typename[k]);
811 /* Add this group as 'SOL' too, for backward compatibility with older gromacs versions */
812 add_grp(gb,gn,nra,aid,"SOL");
816 /* Solvent, create a negated group too */
817 aid=mk_aid(atoms,restype,"Water",&nra,FALSE);
818 if ((nra > 0) && (nra < atoms->nr))
820 add_grp(gb,gn,nra,aid,"non-Water");
827 add_grp(gb,gn,nra,aid,p_typename[k]);
829 analyse_other(restype,atoms,gb,gn,bASK,bVerb);
835 gmx_residuetype_destroy(rt);
837 /* Create a merged water_and_ions group */
843 for(i=0;i<gb->nr;i++)
845 if(!gmx_strcasecmp((*gn)[i],"Water"))
848 nwater = gb->index[i+1]-gb->index[i];
850 else if(!gmx_strcasecmp((*gn)[i],"Ion"))
853 nion = gb->index[i+1]-gb->index[i];
857 if(nwater>0 && nion>0)
859 srenew(gb->index,gb->nr+2);
860 srenew(*gn,gb->nr+1);
861 (*gn)[gb->nr] = strdup("Water_and_ions");
862 srenew(gb->a,gb->nra+nwater+nion);
865 for(i=gb->index[iwater];i<gb->index[iwater+1];i++)
867 gb->a[gb->nra++] = gb->a[i];
872 for(i=gb->index[iion];i<gb->index[iion+1];i++)
874 gb->a[gb->nra++] = gb->a[i];
878 gb->index[gb->nr]=gb->nra;
883 void check_index(char *gname,int n,atom_id index[],char *traj,int natoms)
888 if (index[i] >= natoms)
889 gmx_fatal(FARGS,"%s atom number (index[%d]=%d) is larger than the number of atoms in %s (%d)",
890 gname ? gname : "Index",i+1, index[i]+1,
891 traj ? traj : "the trajectory",natoms);
892 else if (index[i] < 0)
893 gmx_fatal(FARGS,"%s atom number (index[%d]=%d) is less than zero",
894 gname ? gname : "Index",i+1, index[i]+1);
897 t_blocka *init_index(const char *gfile, char ***grpname)
903 char line[STRLEN],*pt,str[STRLEN];
905 in=gmx_fio_fopen(gfile,"r");
907 get_a_line(in,line,STRLEN);
908 if ( line[0]=='[' ) {
917 if (get_header(line,str)) {
919 srenew(b->index,b->nr+1);
920 srenew(*grpname,b->nr);
923 b->index[b->nr]=b->index[b->nr-1];
924 (*grpname)[b->nr-1]=strdup(str);
928 gmx_fatal(FARGS,"The first header of your indexfile is invalid");
931 while (sscanf(pt,"%s",str) == 1) {
935 srenew(b->a,maxentries);
937 b->a[i]=strtol(str, NULL, 10)-1;
940 pt=strstr(pt,str)+strlen(str);
943 } while (get_a_line(in,line,STRLEN));
947 sscanf(line,"%d%d",&b->nr,&b->nra);
948 snew(b->index,b->nr+1);
949 snew(*grpname,b->nr);
952 for (i=0; (i<b->nr); i++) {
953 nread=fscanf(in,"%s%d",str,&ng);
954 (*grpname)[i]=strdup(str);
955 b->index[i+1]=b->index[i]+ng;
956 if (b->index[i+1] > b->nra)
957 gmx_fatal(FARGS,"Something wrong in your indexfile at group %s",str);
958 for(j=0; (j<ng); j++) {
959 nread=fscanf(in,"%d",&a);
960 b->a[b->index[i]+j]=a;
966 for(i=0; (i<b->nr); i++) {
967 for(j=b->index[i]; (j<b->index[i+1]); j++) {
969 fprintf(stderr,"\nWARNING: negative index %d in group %s\n\n",
970 b->a[j],(*grpname)[i]);
977 static void minstring(char *str)
981 for (i=0; (i < (int)strlen(str)); i++)
986 int find_group(char s[], int ngrps, char **grpname)
995 /* first look for whole name match */
997 for(i=0; i<ngrps; i++)
998 if (gmx_strcasecmp_min(s,grpname[i])==0) {
1003 /* second look for first string match */
1005 for(i=0; i<ngrps; i++)
1006 if (gmx_strncasecmp_min(s,grpname[i],n)==0) {
1011 /* last look for arbitrary substring match */
1015 for(i=0; i<ngrps; i++) {
1016 strcpy(string, grpname[i]);
1019 if (strstr(string,s)!=NULL) {
1027 printf("Error: Multiple groups '%s' selected\n", s);
1033 static int qgroup(int *a, int ngrps, char **grpname)
1041 fprintf(stderr,"Select a group: ");
1043 if ( scanf("%s",s)!=1 )
1044 gmx_fatal(FARGS,"Cannot read from input");
1045 trim(s); /* remove spaces */
1046 } while (strlen(s)==0);
1047 aa = strtol(s, &end, 10);
1048 if (aa==0 && end[0] != '\0') /* string entered */
1049 aa = find_group(s, ngrps, grpname);
1050 bInRange = (aa >= 0 && aa < ngrps);
1052 printf("Error: No such group '%s'\n", s);
1053 } while (!bInRange);
1054 printf("Selected %d: '%s'\n", aa, grpname[aa]);
1059 static void rd_groups(t_blocka *grps,char **grpname,char *gnames[],
1060 int ngrps,int isize[],atom_id *index[],int grpnr[])
1065 gmx_fatal(FARGS,"Error: no groups in indexfile");
1066 for(i=0; (i<grps->nr); i++)
1067 fprintf(stderr,"Group %5d (%15s) has %5d elements\n",i,grpname[i],
1068 grps->index[i+1]-grps->index[i]);
1069 for(i=0; (i<ngrps); i++) {
1072 gnr1=qgroup(&grpnr[i], grps->nr, grpname);
1073 if ((gnr1<0) || (gnr1>=grps->nr))
1074 fprintf(stderr,"Select between %d and %d.\n",0,grps->nr-1);
1075 } while ((gnr1<0) || (gnr1>=grps->nr));
1077 fprintf(stderr,"There is one group in the index\n");
1080 gnames[i]=strdup(grpname[gnr1]);
1081 isize[i]=grps->index[gnr1+1]-grps->index[gnr1];
1082 snew(index[i],isize[i]);
1083 for(j=0; (j<isize[i]); j++)
1084 index[i][j]=grps->a[grps->index[gnr1]+j];
1088 void rd_index(const char *statfile,int ngrps,int isize[],
1089 atom_id *index[],char *grpnames[])
1097 gmx_fatal(FARGS,"No index file specified");
1098 grps=init_index(statfile,&gnames);
1099 rd_groups(grps,gnames,grpnames,ngrps,isize,index,grpnr);
1102 void rd_index_nrs(char *statfile,int ngrps,int isize[],
1103 atom_id *index[],char *grpnames[],int grpnr[])
1109 gmx_fatal(FARGS,"No index file specified");
1110 grps=init_index(statfile,&gnames);
1112 rd_groups(grps,gnames,grpnames,ngrps,isize,index,grpnr);
1115 void get_index(t_atoms *atoms, const char *fnm, int ngrps,
1116 int isize[], atom_id *index[],char *grpnames[])
1119 t_blocka *grps = NULL;
1125 grps=init_index(fnm,gnames);
1129 snew(grps->index,1);
1130 analyse(atoms,grps,gnames,FALSE,FALSE);
1133 gmx_incons("You need to supply a valid atoms structure or a valid index file name");
1135 rd_groups(grps,*gnames,grpnames,ngrps,isize,index,grpnr);
1138 t_cluster_ndx *cluster_index(FILE *fplog,const char *ndx)
1144 c->clust = init_index(ndx,&c->grpname);
1146 for(i=0; (i<c->clust->nra); i++)
1147 c->maxframe = max(c->maxframe,c->clust->a[i]);
1148 fprintf(fplog ? fplog : stdout,
1149 "There are %d clusters containing %d structures, highest framenr is %d\n",
1150 c->clust->nr,c->clust->nra,c->maxframe);
1152 pr_blocka(debug,0,"clust",c->clust,TRUE);
1153 for(i=0; (i<c->clust->nra); i++)
1154 if ((c->clust->a[i] < 0) || (c->clust->a[i] > c->maxframe))
1155 gmx_fatal(FARGS,"Range check error for c->clust->a[%d] = %d\n"
1156 "should be within 0 and %d",i,c->clust->a[i],c->maxframe+1);
1158 c->inv_clust=make_invblocka(c->clust,c->maxframe);