1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
53 #include "gmx_fatal.h"
72 #include "fflibutil.h"
84 static const char *res2bb_notermini(const char *name
,
85 int nrr
,const rtprename_t
*rr
)
87 /* NOTE: This function returns the main building block name,
88 * it does not take terminal renaming into account.
93 while (i
< nrr
&& gmx_strcasecmp(name
,rr
[i
].gmx
) != 0) {
97 return (i
< nrr
? rr
[i
].main
: name
);
100 static const char *select_res(int nr
,int resnr
,
101 const char *name
[],const char *expl
[],
103 int nrr
,const rtprename_t
*rr
)
107 printf("Which %s type do you want for residue %d\n",title
,resnr
+1);
108 for(sel
=0; (sel
< nr
); sel
++) {
109 printf("%d. %s (%s)\n",
110 sel
,expl
[sel
],res2bb_notermini(name
[sel
],nrr
,rr
));
112 printf("\nType a number:"); fflush(stdout
);
114 if (scanf("%d",&sel
) != 1)
115 gmx_fatal(FARGS
,"Answer me for res %s %d!",title
,resnr
+1);
120 static const char *get_asptp(int resnr
,int nrr
,const rtprename_t
*rr
)
122 enum { easp
, easpH
, easpNR
};
123 const char *lh
[easpNR
] = { "ASP", "ASPH" };
124 const char *expl
[easpNR
] = {
125 "Not protonated (charge -1)",
126 "Protonated (charge 0)"
129 return select_res(easpNR
,resnr
,lh
,expl
,"ASPARTIC ACID",nrr
,rr
);
132 static const char *get_glutp(int resnr
,int nrr
,const rtprename_t
*rr
)
134 enum { eglu
, egluH
, egluNR
};
135 const char *lh
[egluNR
] = { "GLU", "GLUH" };
136 const char *expl
[egluNR
] = {
137 "Not protonated (charge -1)",
138 "Protonated (charge 0)"
141 return select_res(egluNR
,resnr
,lh
,expl
,"GLUTAMIC ACID",nrr
,rr
);
144 static const char *get_glntp(int resnr
,int nrr
,const rtprename_t
*rr
)
146 enum { egln
, eglnH
, eglnNR
};
147 const char *lh
[eglnNR
] = { "GLN", "QLN" };
148 const char *expl
[eglnNR
] = {
149 "Not protonated (charge 0)",
150 "Protonated (charge +1)"
153 return select_res(eglnNR
,resnr
,lh
,expl
,"GLUTAMINE",nrr
,rr
);
156 static const char *get_lystp(int resnr
,int nrr
,const rtprename_t
*rr
)
158 enum { elys
, elysH
, elysNR
};
159 const char *lh
[elysNR
] = { "LYSN", "LYS" };
160 const char *expl
[elysNR
] = {
161 "Not protonated (charge 0)",
162 "Protonated (charge +1)"
165 return select_res(elysNR
,resnr
,lh
,expl
,"LYSINE",nrr
,rr
);
168 static const char *get_argtp(int resnr
,int nrr
,const rtprename_t
*rr
)
170 enum { earg
, eargH
, eargNR
};
171 const char *lh
[eargNR
] = { "ARGN", "ARG" };
172 const char *expl
[eargNR
] = {
173 "Not protonated (charge 0)",
174 "Protonated (charge +1)"
177 return select_res(eargNR
,resnr
,lh
,expl
,"ARGININE",nrr
,rr
);
180 static const char *get_cystp(int resnr
,int nrr
,const rtprename_t
*rr
)
182 enum { ecys
, ecysH
, ecysNR
};
183 const char *lh
[ecysNR
] = { "CYS2", "CYS" };
184 const char *expl
[ecysNR
] = {
185 "Cysteine in disulfide bridge",
189 return select_res(ecysNR
,resnr
,lh
,expl
,"CYSTEINE",nrr
,rr
);
193 static const char *get_histp(int resnr
,int nrr
,const rtprename_t
*rr
)
195 const char *expl
[ehisNR
] = {
202 return select_res(ehisNR
,resnr
,hh
,expl
,"HISTIDINE",nrr
,rr
);
205 static void read_rtprename(const char *fname
,FILE *fp
,
206 int *nrtprename
,rtprename_t
**rtprename
)
208 char line
[STRLEN
],buf
[STRLEN
];
217 while(get_a_line(fp
,line
,STRLEN
)) {
219 nc
= sscanf(line
,"%s %s %s %s %s %s",
220 rr
[n
].gmx
,rr
[n
].main
,rr
[n
].nter
,rr
[n
].cter
,rr
[n
].bter
,buf
);
222 if (nc
!= 2 && nc
!= 5) {
223 gmx_fatal(FARGS
,"Residue renaming database '%s' has %d columns instead of %d, %d or %d",fname
,ncol
,2,5);
226 } else if (nc
!= ncol
) {
227 gmx_fatal(FARGS
,"A line in residue renaming database '%s' has %d columns, while previous lines have %d columns",fname
,nc
,ncol
);
231 /* This file does not have special termini names, copy them from main */
232 strcpy(rr
[n
].nter
,rr
[n
].main
);
233 strcpy(rr
[n
].cter
,rr
[n
].main
);
234 strcpy(rr
[n
].bter
,rr
[n
].main
);
244 static char *search_resrename(int nrr
,rtprename_t
*rr
,
246 gmx_bool bStart
,gmx_bool bEnd
,
247 gmx_bool bCompareFFRTPname
)
255 while (i
<nrr
&& ((!bCompareFFRTPname
&& strcmp(name
,rr
[i
].gmx
) != 0) ||
256 ( bCompareFFRTPname
&& strcmp(name
,rr
[i
].main
) != 0)))
261 /* If found in the database, rename this residue's rtp building block,
262 * otherwise keep the old name.
284 gmx_fatal(FARGS
,"In the chosen force field there is no residue type for '%s'%s",name
,bStart
? " as a starting terminus" : (bEnd
? " as an ending terminus" : ""));
292 static void rename_resrtp(t_atoms
*pdba
,int nterpairs
,int *r_start
,int *r_end
,
293 int nrr
,rtprename_t
*rr
,t_symtab
*symtab
,
297 gmx_bool bStart
,bEnd
;
299 gmx_bool bFFRTPTERRNM
;
301 bFFRTPTERRNM
= (getenv("GMX_NO_FFRTP_TER_RENAME") == NULL
);
303 for(r
=0; r
<pdba
->nres
; r
++)
307 for(j
=0; j
<nterpairs
; j
++)
314 for(j
=0; j
<nterpairs
; j
++)
322 nn
= search_resrename(nrr
,rr
,*pdba
->resinfo
[r
].rtp
,bStart
,bEnd
,FALSE
);
324 if (bFFRTPTERRNM
&& nn
== NULL
&& (bStart
|| bEnd
))
326 /* This is a terminal residue, but the residue name,
327 * currently stored in .rtp, is not a standard residue name,
328 * but probably a force field specific rtp name.
329 * Check if we need to rename it because it is terminal.
331 nn
= search_resrename(nrr
,rr
,
332 *pdba
->resinfo
[r
].rtp
,bStart
,bEnd
,TRUE
);
335 if (nn
!= NULL
&& strcmp(*pdba
->resinfo
[r
].rtp
,nn
) != 0)
339 printf("Changing rtp entry of residue %d %s to '%s'\n",
340 pdba
->resinfo
[r
].nr
,*pdba
->resinfo
[r
].name
,nn
);
342 pdba
->resinfo
[r
].rtp
= put_symtab(symtab
,nn
);
347 static void pdbres_to_gmxrtp(t_atoms
*pdba
)
351 for(i
=0; (i
<pdba
->nres
); i
++)
353 if (pdba
->resinfo
[i
].rtp
== NULL
)
355 pdba
->resinfo
[i
].rtp
= pdba
->resinfo
[i
].name
;
360 static void rename_pdbres(t_atoms
*pdba
,const char *oldnm
,const char *newnm
,
361 gmx_bool bFullCompare
,t_symtab
*symtab
)
366 for(i
=0; (i
<pdba
->nres
); i
++)
368 resnm
= *pdba
->resinfo
[i
].name
;
369 if ((bFullCompare
&& (gmx_strcasecmp(resnm
,oldnm
) == 0)) ||
370 (!bFullCompare
&& strstr(resnm
,oldnm
) != NULL
))
372 /* Rename the residue name (not the rtp name) */
373 pdba
->resinfo
[i
].name
= put_symtab(symtab
,newnm
);
378 static void rename_bb(t_atoms
*pdba
,const char *oldnm
,const char *newnm
,
379 gmx_bool bFullCompare
,t_symtab
*symtab
)
384 for(i
=0; (i
<pdba
->nres
); i
++)
386 /* We have not set the rtp name yes, use the residue name */
387 bbnm
= *pdba
->resinfo
[i
].name
;
388 if ((bFullCompare
&& (gmx_strcasecmp(bbnm
,oldnm
) == 0)) ||
389 (!bFullCompare
&& strstr(bbnm
,oldnm
) != NULL
))
391 /* Change the rtp builing block name */
392 pdba
->resinfo
[i
].rtp
= put_symtab(symtab
,newnm
);
397 static void rename_bbint(t_atoms
*pdba
,const char *oldnm
,
398 const char *gettp(int,int,const rtprename_t
*),
399 gmx_bool bFullCompare
,
401 int nrr
,const rtprename_t
*rr
)
407 for(i
=0; i
<pdba
->nres
; i
++)
409 /* We have not set the rtp name yes, use the residue name */
410 bbnm
= *pdba
->resinfo
[i
].name
;
411 if ((bFullCompare
&& (strcmp(bbnm
,oldnm
) == 0)) ||
412 (!bFullCompare
&& strstr(bbnm
,oldnm
) != NULL
))
414 ptr
= gettp(i
,nrr
,rr
);
415 pdba
->resinfo
[i
].rtp
= put_symtab(symtab
,ptr
);
420 static void check_occupancy(t_atoms
*atoms
,const char *filename
,gmx_bool bVerbose
)
426 ftp
= fn2ftp(filename
);
427 if (!atoms
->pdbinfo
|| ((ftp
!= efPDB
) && (ftp
!= efBRK
) && (ftp
!= efENT
)))
428 fprintf(stderr
,"No occupancies in %s\n",filename
);
430 for(i
=0; (i
<atoms
->nr
); i
++) {
431 if (atoms
->pdbinfo
[i
].occup
!= 1) {
433 fprintf(stderr
,"Occupancy for atom %s%d-%s is %f rather than 1\n",
434 *atoms
->resinfo
[atoms
->atom
[i
].resind
].name
,
435 atoms
->resinfo
[ atoms
->atom
[i
].resind
].nr
,
437 atoms
->pdbinfo
[i
].occup
);
438 if (atoms
->pdbinfo
[i
].occup
== 0)
444 if (nzero
== atoms
->nr
) {
445 fprintf(stderr
,"All occupancy fields zero. This is probably not an X-Ray structure\n");
446 } else if ((nzero
> 0) || (nnotone
> 0)) {
449 "WARNING: there were %d atoms with zero occupancy and %d atoms with\n"
450 " occupancy unequal to one (out of %d atoms). Check your pdb file.\n"
452 nzero
,nnotone
,atoms
->nr
);
454 fprintf(stderr
,"All occupancies are one\n");
459 void write_posres(char *fn
,t_atoms
*pdba
,real fc
)
464 fp
=gmx_fio_fopen(fn
,"w");
466 "; In this topology include file, you will find position restraint\n"
467 "; entries for all the heavy atoms in your original pdb file.\n"
468 "; This means that all the protons which were added by pdb2gmx are\n"
469 "; not restrained.\n"
471 "[ position_restraints ]\n"
472 "; %4s%6s%8s%8s%8s\n","atom","type","fx","fy","fz"
474 for(i
=0; (i
<pdba
->nr
); i
++) {
475 if (!is_hydrogen(*pdba
->atomname
[i
]) && !is_dummymass(*pdba
->atomname
[i
]))
476 fprintf(fp
,"%6d%6d %g %g %g\n",i
+1,1,fc
,fc
,fc
);
481 static int read_pdball(const char *inf
, const char *outf
,char *title
,
482 t_atoms
*atoms
, rvec
**x
,
483 int *ePBC
,matrix box
, gmx_bool bRemoveH
,
484 t_symtab
*symtab
,gmx_residuetype_t rt
,const char *watres
,
485 gmx_atomprop_t aps
,gmx_bool bVerbose
)
486 /* Read a pdb file. (containing proteins) */
488 int natom
,new_natom
,i
;
491 printf("Reading %s...\n",inf
);
492 get_stx_coordnum(inf
,&natom
);
493 init_t_atoms(atoms
,natom
,TRUE
);
495 read_stx_conf(inf
,title
,atoms
,*x
,NULL
,ePBC
,box
);
496 if (fn2ftp(inf
) == efPDB
)
497 get_pdb_atomnumber(atoms
,aps
);
500 for(i
=0; i
<atoms
->nr
; i
++)
501 if (!is_hydrogen(*atoms
->atomname
[i
])) {
502 atoms
->atom
[new_natom
]=atoms
->atom
[i
];
503 atoms
->atomname
[new_natom
]=atoms
->atomname
[i
];
504 atoms
->pdbinfo
[new_natom
]=atoms
->pdbinfo
[i
];
505 copy_rvec((*x
)[i
],(*x
)[new_natom
]);
513 if (title
&& title
[0])
514 printf(" '%s',",title
);
515 printf(" %d atoms\n",natom
);
517 /* Rename residues */
518 rename_pdbres(atoms
,"HOH",watres
,FALSE
,symtab
);
519 rename_pdbres(atoms
,"SOL",watres
,FALSE
,symtab
);
520 rename_pdbres(atoms
,"WAT",watres
,FALSE
,symtab
);
522 rename_atoms("xlateat.dat",NULL
,
523 atoms
,symtab
,NULL
,TRUE
,rt
,TRUE
,bVerbose
);
529 write_sto_conf(outf
,title
,atoms
,*x
,NULL
,*ePBC
,box
);
534 void process_chain(t_atoms
*pdba
, rvec
*x
,
535 gmx_bool bTrpU
,gmx_bool bPheU
,gmx_bool bTyrU
,
536 gmx_bool bLysMan
,gmx_bool bAspMan
,gmx_bool bGluMan
,
537 gmx_bool bHisMan
,gmx_bool bArgMan
,gmx_bool bGlnMan
,
538 real angle
,real distance
,t_symtab
*symtab
,
539 int nrr
,const rtprename_t
*rr
)
541 /* Rename aromatics, lys, asp and histidine */
542 if (bTyrU
) rename_bb(pdba
,"TYR","TYRU",FALSE
,symtab
);
543 if (bTrpU
) rename_bb(pdba
,"TRP","TRPU",FALSE
,symtab
);
544 if (bPheU
) rename_bb(pdba
,"PHE","PHEU",FALSE
,symtab
);
546 rename_bbint(pdba
,"LYS",get_lystp
,FALSE
,symtab
,nrr
,rr
);
548 rename_bbint(pdba
,"ARG",get_argtp
,FALSE
,symtab
,nrr
,rr
);
550 rename_bbint(pdba
,"GLN",get_glntp
,FALSE
,symtab
,nrr
,rr
);
552 rename_bbint(pdba
,"ASP",get_asptp
,FALSE
,symtab
,nrr
,rr
);
554 rename_bb(pdba
,"ASPH","ASP",FALSE
,symtab
);
556 rename_bbint(pdba
,"GLU",get_glutp
,FALSE
,symtab
,nrr
,rr
);
558 rename_bb(pdba
,"GLUH","GLU",FALSE
,symtab
);
561 set_histp(pdba
,x
,angle
,distance
);
563 rename_bbint(pdba
,"HIS",get_histp
,TRUE
,symtab
,nrr
,rr
);
565 /* Initialize the rtp builing block names with the residue names
566 * for the residues that have not been processed above.
568 pdbres_to_gmxrtp(pdba
);
570 /* Now we have all rtp names set.
571 * The rtp names will conform to Gromacs naming,
572 * unless the input pdb file contained one or more force field specific
573 * rtp names as residue names.
577 /* struct for sorting the atoms from the pdb file */
579 int resnr
; /* residue number */
580 int j
; /* database order index */
581 int index
; /* original atom number */
582 char anm1
; /* second letter of atom name */
583 char altloc
; /* alternate location indicator */
586 int pdbicomp(const void *a
,const void *b
)
594 d
= (pa
->resnr
- pb
->resnr
);
598 d
= (pa
->anm1
- pb
->anm1
);
600 d
= (pa
->altloc
- pb
->altloc
);
607 static void sort_pdbatoms(int nrtp
,t_restp restp
[],t_hackblock hb
[],
608 int natoms
,t_atoms
**pdbaptr
,rvec
**x
,
609 t_blocka
*block
,char ***gnames
)
611 t_atoms
*pdba
,*pdbnew
;
626 for(i
=0; i
<natoms
; i
++)
628 atomnm
= *pdba
->atomname
[i
];
629 rptr
= &restp
[pdba
->atom
[i
].resind
];
630 for(j
=0; (j
<rptr
->natom
); j
++)
632 if (gmx_strcasecmp(atomnm
,*(rptr
->atomname
[j
])) == 0)
642 "Atom %s in residue %s %d was not found in rtp entry %s with %d atoms\n"
643 "while sorting atoms.\n%s",atomnm
,
644 *pdba
->resinfo
[pdba
->atom
[i
].resind
].name
,
645 pdba
->resinfo
[pdba
->atom
[i
].resind
].nr
,
648 is_hydrogen(atomnm
) ?
649 "\nFor a hydrogen, this can be a different protonation state, or it\n"
650 "might have had a different number in the PDB file and was rebuilt\n"
651 "(it might for instance have been H3, and we only expected H1 & H2).\n"
652 "Note that hydrogens might have been added to the entry for the N-terminus.\n"
653 "Remove this hydrogen or choose a different protonation state to solve it.\n"
654 "Option -ignh will ignore all hydrogens in the input." : ".");
655 gmx_fatal(FARGS
,buf
);
657 /* make shadow array to be sorted into indexgroup */
658 pdbi
[i
].resnr
= pdba
->atom
[i
].resind
;
661 pdbi
[i
].anm1
= atomnm
[1];
662 pdbi
[i
].altloc
= pdba
->pdbinfo
[i
].altloc
;
664 qsort(pdbi
,natoms
,(size_t)sizeof(pdbi
[0]),pdbicomp
);
666 /* pdba is sorted in pdbnew using the pdbi index */
669 init_t_atoms(pdbnew
,natoms
,TRUE
);
672 pdbnew
->nres
=pdba
->nres
;
673 sfree(pdbnew
->resinfo
);
674 pdbnew
->resinfo
=pdba
->resinfo
;
675 for (i
=0; i
<natoms
; i
++) {
676 pdbnew
->atom
[i
] = pdba
->atom
[pdbi
[i
].index
];
677 pdbnew
->atomname
[i
] = pdba
->atomname
[pdbi
[i
].index
];
678 pdbnew
->pdbinfo
[i
] = pdba
->pdbinfo
[pdbi
[i
].index
];
679 copy_rvec((*x
)[pdbi
[i
].index
],(*xnew
)[i
]);
680 /* make indexgroup in block */
684 sfree(pdba
->atomname
);
686 sfree(pdba
->pdbinfo
);
689 /* copy the sorted pdbnew back to pdba */
692 add_grp(block
, gnames
, natoms
, a
, "prot_sort");
698 static int remove_duplicate_atoms(t_atoms
*pdba
,rvec x
[],gmx_bool bVerbose
)
700 int i
,j
,oldnatoms
,ndel
;
703 printf("Checking for duplicate atoms....\n");
704 oldnatoms
= pdba
->nr
;
706 /* NOTE: pdba->nr is modified inside the loop */
707 for(i
=1; (i
< pdba
->nr
); i
++) {
708 /* compare 'i' and 'i-1', throw away 'i' if they are identical
709 this is a 'while' because multiple alternate locations can be present */
710 while ( (i
< pdba
->nr
) &&
711 (pdba
->atom
[i
-1].resind
== pdba
->atom
[i
].resind
) &&
712 (strcmp(*pdba
->atomname
[i
-1],*pdba
->atomname
[i
])==0) ) {
715 ri
= &pdba
->resinfo
[pdba
->atom
[i
].resind
];
716 printf("deleting duplicate atom %4s %s%4d%c",
717 *pdba
->atomname
[i
],*ri
->name
,ri
->nr
,ri
->ic
);
718 if (ri
->chainid
&& (ri
->chainid
!= ' '))
719 printf(" ch %c", ri
->chainid
);
721 if (pdba
->pdbinfo
[i
].atomnr
)
722 printf(" pdb nr %4d",pdba
->pdbinfo
[i
].atomnr
);
723 if (pdba
->pdbinfo
[i
].altloc
&& (pdba
->pdbinfo
[i
].altloc
!=' '))
724 printf(" altloc %c",pdba
->pdbinfo
[i
].altloc
);
729 /* We can not free, since it might be in the symtab */
730 /* sfree(pdba->atomname[i]); */
731 for (j
=i
; j
< pdba
->nr
; j
++) {
732 pdba
->atom
[j
] = pdba
->atom
[j
+1];
733 pdba
->atomname
[j
] = pdba
->atomname
[j
+1];
734 pdba
->pdbinfo
[j
] = pdba
->pdbinfo
[j
+1];
735 copy_rvec(x
[j
+1],x
[j
]);
737 srenew(pdba
->atom
, pdba
->nr
);
738 /* srenew(pdba->atomname, pdba->nr); */
739 srenew(pdba
->pdbinfo
, pdba
->nr
);
742 if (pdba
->nr
!= oldnatoms
)
743 printf("Now there are %d atoms. Deleted %d duplicates.\n",pdba
->nr
,ndel
);
748 void find_nc_ter(t_atoms
*pdba
,int r0
,int r1
,int *r_start
,int *r_end
,gmx_residuetype_t rt
)
751 const char *p_startrestype
;
752 const char *p_restype
;
753 int nstartwarn
,nendwarn
;
761 /* Find the starting terminus (typially N or 5') */
762 for(i
=r0
;i
<r1
&& *r_start
==-1;i
++)
764 gmx_residuetype_get_type(rt
,*pdba
->resinfo
[i
].name
,&p_startrestype
);
765 if( !gmx_strcasecmp(p_startrestype
,"Protein") || !gmx_strcasecmp(p_startrestype
,"DNA") || !gmx_strcasecmp(p_startrestype
,"RNA") )
767 printf("Identified residue %s%d as a starting terminus.\n",*pdba
->resinfo
[i
].name
,pdba
->resinfo
[i
].nr
);
774 printf("Warning: Starting residue %s%d in chain not identified as Protein/RNA/DNA.\n",*pdba
->resinfo
[i
].name
,pdba
->resinfo
[i
].nr
);
778 printf("More than 5 unidentified residues at start of chain - disabling further warnings.\n");
786 /* Go through the rest of the residues, check that they are the same class, and identify the ending terminus. */
787 for(i
=*r_start
;i
<r1
;i
++)
789 gmx_residuetype_get_type(rt
,*pdba
->resinfo
[i
].name
,&p_restype
);
790 if( !gmx_strcasecmp(p_restype
,p_startrestype
) && nendwarn
==0)
798 printf("Warning: Residue %s%d in chain has different type (%s) from starting residue %s%d (%s).\n",
799 *pdba
->resinfo
[i
].name
,pdba
->resinfo
[i
].nr
,p_restype
,
800 *pdba
->resinfo
[*r_start
].name
,pdba
->resinfo
[*r_start
].nr
,p_startrestype
);
804 printf("More than 5 unidentified residues at end of chain - disabling further warnings.\n");
813 printf("Identified residue %s%d as a ending terminus.\n",*pdba
->resinfo
[*r_end
].name
,pdba
->resinfo
[*r_end
].nr
);
819 modify_chain_numbers(t_atoms
* pdba
,
820 const char * chainsep
)
823 char old_prev_chainid
;
824 char old_this_chainid
;
825 int old_prev_chainnum
;
826 int old_this_chainnum
;
839 printf("Splitting PDB chains based on ");
840 splitting
= SPLIT_TER_ONLY
; /* keep compiler happy */
842 /* Be a bit flexible to catch typos */
843 if (!strncmp(chainsep
,"id_o",4) || !strncmp(chainsep
,"int",3))
845 /* For later interactive splitting we tentatively assign new chain numbers at either changing id or ter records */
846 splitting
= SPLIT_ID_OR_TER
;
847 printf("TER records or changing chain id.\n");
849 else if (!strncmp(chainsep
,"id_a",4))
851 splitting
= SPLIT_ID_AND_TER
;
852 printf("TER records and chain id.\n");
854 else if (strlen(chainsep
)==2 && !strncmp(chainsep
,"id",4))
856 splitting
= SPLIT_ID_ONLY
;
857 printf("changing chain id only (ignoring TER records).\n");
859 else if (chainsep
[0]=='t')
861 splitting
= SPLIT_TER_ONLY
;
862 printf("TER records only (ignoring chain id).\n");
866 gmx_fatal(FARGS
,"Unidentified setting for chain separation: %s\n",chainsep
);
869 /* The default chain enumeration is based on TER records only, which is reflected in chainnum below */
871 old_prev_chainid
= '?';
872 old_prev_chainnum
= -1;
875 for(i
=0;i
<pdba
->nres
;i
++)
877 ri
= &pdba
->resinfo
[i
];
878 old_this_chainid
= ri
->chainid
;
879 old_this_chainnum
= ri
->chainnum
;
883 case SPLIT_ID_OR_TER
:
884 if(old_this_chainid
!= old_prev_chainid
|| old_this_chainnum
!= old_prev_chainnum
)
890 case SPLIT_ID_AND_TER
:
891 if(old_this_chainid
!= old_prev_chainid
&& old_this_chainnum
!= old_prev_chainnum
)
898 if(old_this_chainid
!= old_prev_chainid
)
905 if(old_this_chainnum
!= old_prev_chainnum
)
911 gmx_fatal(FARGS
,"Internal inconsistency - this shouldn't happen...");
914 old_prev_chainid
= old_this_chainid
;
915 old_prev_chainnum
= old_this_chainnum
;
917 ri
->chainnum
= new_chainnum
;
946 int main(int argc
, char *argv
[])
948 const char *desc
[] = {
949 "This program reads a pdb (or gro) file, reads",
950 "some database files, adds hydrogens to the molecules and generates",
951 "coordinates in Gromacs (Gromos), or optionally pdb, format",
952 "and a topology in Gromacs format.",
953 "These files can subsequently be processed to generate a run input file.",
955 "pdb2gmx will search for force fields by looking for",
956 "a [TT]forcefield.itp[tt] file in subdirectories [TT]<forcefield>.ff[tt]",
957 "of the current working directory and of the Gromacs library directory",
958 "as inferred from the path of the binary or the [TT]GMXLIB[tt] environment",
960 "By default the forcefield selection is interactive,",
961 "but you can use the [TT]-ff[tt] option to specify one of the short names",
962 "in the list on the command line instead. In that case pdb2gmx just looks",
963 "for the corresponding [TT]<forcefield>.ff[tt] directory.",
965 "After choosing a force field, all files will be read only from",
966 "the corresponding force field directory.",
967 "If you want to modify or add a residue types, you can copy the force",
968 "field directory from the Gromacs library directory to your current",
969 "working directory. If you want to add new protein residue types,",
970 "you will need to modify residuetypes.dat in the library directory",
971 "or copy the whole library directory to a local directory and set",
972 "the environment variable [TT]GMXLIB[tt] to the name of that directory.",
973 "Check chapter 5 of the manual for more information about file formats.",
976 "Note that a pdb file is nothing more than a file format, and it",
977 "need not necessarily contain a protein structure. Every kind of",
978 "molecule for which there is support in the database can be converted.",
979 "If there is no support in the database, you can add it yourself.[PAR]",
981 "The program has limited intelligence, it reads a number of database",
982 "files, that allow it to make special bonds (Cys-Cys, Heme-His, etc.),",
983 "if necessary this can be done manually. The program can prompt the",
984 "user to select which kind of LYS, ASP, GLU, CYS or HIS residue she",
985 "wants. For LYS the choice is between neutral (two protons on NZ) or",
986 "protonated (three protons, default), for ASP and GLU unprotonated",
987 "(default) or protonated, for HIS the proton can be either on ND1,",
988 "on NE2 or on both. By default these selections are done automatically.",
989 "For His, this is based on an optimal hydrogen bonding",
990 "conformation. Hydrogen bonds are defined based on a simple geometric",
991 "criterion, specified by the maximum hydrogen-donor-acceptor angle",
992 "and donor-acceptor distance, which are set by [TT]-angle[tt] and",
993 "[TT]-dist[tt] respectively.[PAR]",
995 "The separation of chains is not entirely trivial since the markup",
996 "in user-generated PDB files frequently varies and sometimes it",
997 "is desirable to merge entries across a TER record, for instance",
998 "if you want a disulfide bridge or distance restraints between",
999 "two protein chains or if you have a HEME group bound to a protein.",
1000 "In such cases multiple chains should be contained in a single",
1001 "[TT]molecule_type[tt] definition.",
1002 "To handle this, pdb2gmx has an option [TT]-chainsep[tt] so you can",
1003 "choose whether a new chain should start when we find a TER record,",
1004 "when the chain id changes, combinations of either or both of these",
1005 "or fully interactively.[PAR]",
1007 "pdb2gmx will also check the occupancy field of the pdb file.",
1008 "If any of the occupancies are not one, indicating that the atom is",
1009 "not resolved well in the structure, a warning message is issued.",
1010 "When a pdb file does not originate from an X-Ray structure determination",
1011 "all occupancy fields may be zero. Either way, it is up to the user",
1012 "to verify the correctness of the input data (read the article!).[PAR]",
1014 "During processing the atoms will be reordered according to Gromacs",
1015 "conventions. With [TT]-n[tt] an index file can be generated that",
1016 "contains one group reordered in the same way. This allows you to",
1017 "convert a Gromos trajectory and coordinate file to Gromos. There is",
1018 "one limitation: reordering is done after the hydrogens are stripped",
1019 "from the input and before new hydrogens are added. This means that",
1020 "you should not use [TT]-ignh[tt].[PAR]",
1022 "The [TT].gro[tt] and [TT].g96[tt] file formats do not support chain",
1023 "identifiers. Therefore it is useful to enter a pdb file name at",
1024 "the [TT]-o[tt] option when you want to convert a multi-chain pdb file.",
1027 "The option [TT]-vsite[tt] removes hydrogen and fast improper dihedral",
1028 "motions. Angular and out-of-plane motions can be removed by changing",
1029 "hydrogens into virtual sites and fixing angles, which fixes their",
1030 "position relative to neighboring atoms. Additionally, all atoms in the",
1031 "aromatic rings of the standard amino acids (i.e. PHE, TRP, TYR and HIS)",
1032 "can be converted into virtual sites, eliminating the fast improper dihedral",
1033 "fluctuations in these rings. Note that in this case all other hydrogen",
1034 "atoms are also converted to virtual sites. The mass of all atoms that are",
1035 "converted into virtual sites, is added to the heavy atoms.[PAR]",
1036 "Also slowing down of dihedral motion can be done with [TT]-heavyh[tt]",
1037 "done by increasing the hydrogen-mass by a factor of 4. This is also",
1038 "done for water hydrogens to slow down the rotational motion of water.",
1039 "The increase in mass of the hydrogens is subtracted from the bonded",
1040 "(heavy) atom so that the total mass of the system remains the same."
1044 FILE *fp
,*top_file
,*top_file2
,*itp_file
=NULL
;
1046 t_atoms pdba_all
,*pdba
;
1050 int chain
,nch
,maxch
,nwaterchain
;
1052 t_chain
*chains
,*cc
;
1053 char select
[STRLEN
];
1066 gpp_atomtype_t atype
;
1067 gmx_residuetype_t rt
;
1069 char fn
[256],itp_fn
[STRLEN
],posre_fn
[STRLEN
],buf_fn
[STRLEN
];
1070 char molname
[STRLEN
],title
[STRLEN
],quote
[STRLEN
],generator
[STRLEN
];
1071 char *c
,forcefield
[STRLEN
],ffdir
[STRLEN
];
1072 char ffname
[STRLEN
],suffix
[STRLEN
],buf
[STRLEN
];
1081 rtprename_t
*rtprename
=NULL
;
1082 int nah
,nNtdb
,nCtdb
,ntdblist
;
1083 t_hackblock
*ntdb
,*ctdb
,**tdblist
;
1087 gmx_bool bVsites
=FALSE
,bWat
,bPrevWat
=FALSE
,bITP
,bVsiteAromatics
=FALSE
,bMerge
;
1089 t_hackblock
*hb_chain
;
1090 t_restp
*restp_chain
;
1092 const char *p_restype
;
1096 const char * prev_atomname
;
1097 const char * this_atomname
;
1098 const char * prev_resname
;
1099 const char * this_resname
;
1104 int prev_chainnumber
;
1105 int this_chainnumber
;
1107 int this_chainstart
;
1108 int prev_chainstart
;
1113 { efSTX
, "-f", "eiwit.pdb", ffREAD
},
1114 { efSTO
, "-o", "conf", ffWRITE
},
1115 { efTOP
, NULL
, NULL
, ffWRITE
},
1116 { efITP
, "-i", "posre", ffWRITE
},
1117 { efNDX
, "-n", "clean", ffOPTWR
},
1118 { efSTO
, "-q", "clean.pdb", ffOPTWR
}
1120 #define NFILE asize(fnm)
1123 /* Command line arguments must be static */
1124 static gmx_bool bNewRTP
=FALSE
;
1125 static gmx_bool bInter
=FALSE
, bCysMan
=FALSE
;
1126 static gmx_bool bLysMan
=FALSE
, bAspMan
=FALSE
, bGluMan
=FALSE
, bHisMan
=FALSE
;
1127 static gmx_bool bGlnMan
=FALSE
, bArgMan
=FALSE
;
1128 static gmx_bool bTerMan
=FALSE
, bUnA
=FALSE
, bHeavyH
;
1129 static gmx_bool bSort
=TRUE
, bAllowMissing
=FALSE
, bRemoveH
=FALSE
;
1130 static gmx_bool bDeuterate
=FALSE
,bVerbose
=FALSE
,bChargeGroups
=TRUE
,bCmap
=TRUE
;
1131 static gmx_bool bRenumRes
=FALSE
,bRTPresname
=FALSE
;
1132 static real angle
=135.0, distance
=0.3,posre_fc
=1000;
1133 static real long_bond_dist
=0.25, short_bond_dist
=0.05;
1134 static const char *vsitestr
[] = { NULL
, "none", "hydrogens", "aromatics", NULL
};
1135 static const char *watstr
[] = { NULL
, "select", "none", "spc", "spce", "tip3p", "tip4p", "tip5p", NULL
};
1136 static const char *chainsep
[] = { NULL
, "id_or_ter", "id_and_ter", "ter", "id", "interactive", NULL
};
1137 static const char *ff
= "select";
1140 { "-newrtp", FALSE
, etBOOL
, {&bNewRTP
},
1141 "HIDDENWrite the residue database in new format to 'new.rtp'"},
1142 { "-lb", FALSE
, etREAL
, {&long_bond_dist
},
1143 "HIDDENLong bond warning distance" },
1144 { "-sb", FALSE
, etREAL
, {&short_bond_dist
},
1145 "HIDDENShort bond warning distance" },
1146 { "-chainsep", FALSE
, etENUM
, {chainsep
},
1147 "Condition in PDB files when a new chain and molecule_type should be started" },
1148 { "-ff", FALSE
, etSTR
, {&ff
},
1149 "Force field, interactive by default. Use -h for information." },
1150 { "-water", FALSE
, etENUM
, {watstr
},
1151 "Water model to use" },
1152 { "-inter", FALSE
, etBOOL
, {&bInter
},
1153 "Set the next 8 options to interactive"},
1154 { "-ss", FALSE
, etBOOL
, {&bCysMan
},
1155 "Interactive SS bridge selection" },
1156 { "-ter", FALSE
, etBOOL
, {&bTerMan
},
1157 "Interactive termini selection, iso charged" },
1158 { "-lys", FALSE
, etBOOL
, {&bLysMan
},
1159 "Interactive Lysine selection, iso charged" },
1160 { "-arg", FALSE
, etBOOL
, {&bArgMan
},
1161 "Interactive Arganine selection, iso charged" },
1162 { "-asp", FALSE
, etBOOL
, {&bAspMan
},
1163 "Interactive Aspartic Acid selection, iso charged" },
1164 { "-glu", FALSE
, etBOOL
, {&bGluMan
},
1165 "Interactive Glutamic Acid selection, iso charged" },
1166 { "-gln", FALSE
, etBOOL
, {&bGlnMan
},
1167 "Interactive Glutamine selection, iso neutral" },
1168 { "-his", FALSE
, etBOOL
, {&bHisMan
},
1169 "Interactive Histidine selection, iso checking H-bonds" },
1170 { "-angle", FALSE
, etREAL
, {&angle
},
1171 "Minimum hydrogen-donor-acceptor angle for a H-bond (degrees)" },
1172 { "-dist", FALSE
, etREAL
, {&distance
},
1173 "Maximum donor-acceptor distance for a H-bond (nm)" },
1174 { "-una", FALSE
, etBOOL
, {&bUnA
},
1175 "Select aromatic rings with united CH atoms on Phenylalanine, "
1176 "Tryptophane and Tyrosine" },
1177 { "-sort", FALSE
, etBOOL
, {&bSort
},
1178 "HIDDENSort the residues according to database, turning this off is dangerous as charge groups might be broken in parts" },
1179 { "-ignh", FALSE
, etBOOL
, {&bRemoveH
},
1180 "Ignore hydrogen atoms that are in the pdb file" },
1181 { "-missing",FALSE
, etBOOL
, {&bAllowMissing
},
1182 "Continue when atoms are missing, dangerous" },
1183 { "-v", FALSE
, etBOOL
, {&bVerbose
},
1184 "Be slightly more verbose in messages" },
1185 { "-posrefc",FALSE
, etREAL
, {&posre_fc
},
1186 "Force constant for position restraints" },
1187 { "-vsite", FALSE
, etENUM
, {vsitestr
},
1188 "Convert atoms to virtual sites" },
1189 { "-heavyh", FALSE
, etBOOL
, {&bHeavyH
},
1190 "Make hydrogen atoms heavy" },
1191 { "-deuterate", FALSE
, etBOOL
, {&bDeuterate
},
1192 "Change the mass of hydrogens to 2 amu" },
1193 { "-chargegrp", TRUE
, etBOOL
, {&bChargeGroups
},
1194 "Use charge groups in the rtp file" },
1195 { "-cmap", TRUE
, etBOOL
, {&bCmap
},
1196 "Use cmap torsions (if enabled in the rtp file)" },
1197 { "-renum", TRUE
, etBOOL
, {&bRenumRes
},
1198 "Renumber the residues consecutively in the output" },
1199 { "-rtpres", TRUE
, etBOOL
, {&bRTPresname
},
1200 "Use rtp entry names as residue names" }
1202 #define NPARGS asize(pa)
1204 CopyRight(stderr
,argv
[0]);
1205 parse_common_args(&argc
,argv
,0,NFILE
,fnm
,asize(pa
),pa
,asize(desc
),desc
,
1208 /* Force field selection, interactive or direct */
1209 choose_ff(strcmp(ff
,"select") == 0 ? NULL
: ff
,
1210 forcefield
,sizeof(forcefield
),
1211 ffdir
,sizeof(ffdir
));
1213 if (strlen(forcefield
) > 0) {
1214 strcpy(ffname
,forcefield
);
1215 ffname
[0] = toupper(ffname
[0]);
1217 gmx_fatal(FARGS
,"Empty forcefield string");
1220 printf("\nUsing the %s force field in directory %s\n\n",
1223 choose_watermodel(watstr
[0],ffdir
,&watermodel
);
1226 /* if anything changes here, also change description of -inter */
1239 else if (bDeuterate
)
1244 switch(vsitestr
[0][0]) {
1245 case 'n': /* none */
1247 bVsiteAromatics
=FALSE
;
1249 case 'h': /* hydrogens */
1251 bVsiteAromatics
=FALSE
;
1253 case 'a': /* aromatics */
1255 bVsiteAromatics
=TRUE
;
1258 gmx_fatal(FARGS
,"DEATH HORROR in $s (%d): vsitestr[0]='%s'",
1259 __FILE__
,__LINE__
,vsitestr
[0]);
1262 /* Open the symbol table */
1263 open_symtab(&symtab
);
1265 /* Residue type database */
1266 gmx_residuetype_init(&rt
);
1268 /* Read residue renaming database(s), if present */
1269 nrrn
= fflib_search_file_end(ffdir
,".r2b",FALSE
,&rrn
);
1273 for(i
=0; i
<nrrn
; i
++) {
1274 fp
= fflib_open(rrn
[i
]);
1275 read_rtprename(rrn
[i
],fp
,&nrtprename
,&rtprename
);
1281 /* Add all alternative names from the residue renaming database to the list of recognized amino/nucleic acids. */
1283 for(i
=0;i
<nrtprename
;i
++)
1285 rc
=gmx_residuetype_get_type(rt
,rtprename
[i
].gmx
,&p_restype
);
1287 /* Only add names if the 'standard' gromacs/iupac base name was found */
1290 gmx_residuetype_add(rt
,rtprename
[i
].main
,p_restype
);
1291 gmx_residuetype_add(rt
,rtprename
[i
].nter
,p_restype
);
1292 gmx_residuetype_add(rt
,rtprename
[i
].cter
,p_restype
);
1293 gmx_residuetype_add(rt
,rtprename
[i
].bter
,p_restype
);
1298 if (watermodel
!= NULL
&& (strstr(watermodel
,"4p") ||
1299 strstr(watermodel
,"4P"))) {
1301 } else if (watermodel
!= NULL
&& (strstr(watermodel
,"5p") ||
1302 strstr(watermodel
,"5P"))) {
1308 aps
= gmx_atomprop_init();
1309 natom
= read_pdball(opt2fn("-f",NFILE
,fnm
),opt2fn_null("-q",NFILE
,fnm
),title
,
1310 &pdba_all
,&pdbx
,&ePBC
,box
,bRemoveH
,&symtab
,rt
,watres
,
1314 gmx_fatal(FARGS
,"No atoms found in pdb file %s\n",opt2fn("-f",NFILE
,fnm
));
1316 printf("Analyzing pdb file\n");
1321 modify_chain_numbers(&pdba_all
,chainsep
[0]);
1324 bMerge
= !strncmp(chainsep
[0],"int",3);
1326 this_atomname
= NULL
;
1328 this_resname
= NULL
;
1331 this_chainnumber
= -1;
1332 this_chainstart
= 0;
1335 for (i
=0; (i
<natom
); i
++)
1337 ri
= &pdba_all
.resinfo
[pdba_all
.atom
[i
].resind
];
1339 prev_atomname
= this_atomname
;
1340 prev_atomnum
= this_atomnum
;
1341 prev_resname
= this_resname
;
1342 prev_resnum
= this_resnum
;
1343 prev_chainid
= this_chainid
;
1344 prev_chainnumber
= this_chainnumber
;
1345 prev_chainstart
= this_chainstart
;
1347 this_atomname
= *pdba_all
.atomname
[i
];
1348 this_atomnum
= (pdba_all
.pdbinfo
!= NULL
) ? pdba_all
.pdbinfo
[i
].atomnr
: i
+1;
1349 this_resname
= *ri
->name
;
1350 this_resnum
= ri
->nr
;
1351 this_chainid
= ri
->chainid
;
1352 this_chainnumber
= ri
->chainnum
;
1354 bWat
= gmx_strcasecmp(*ri
->name
,watres
) == 0;
1355 if ((i
== 0) || (this_chainnumber
!= prev_chainnumber
) || (bWat
!= bPrevWat
))
1357 this_chainstart
= pdba_all
.atom
[i
].resind
;
1358 if (bMerge
&& i
>0 && !bWat
)
1360 printf("Merge chain ending with residue %s%d (chain id '%c', atom %d %s) with\n"
1361 "chain starting with residue %s%d (chain id '%c', atom %d %s)? [n/y]\n",
1362 prev_resname
,prev_resnum
,prev_chainid
,prev_atomnum
,prev_atomname
,
1363 this_resname
,this_resnum
,this_chainid
,this_atomnum
,this_atomname
);
1365 if(NULL
==fgets(select
,STRLEN
-1,stdin
))
1367 gmx_fatal(FARGS
,"Error reading from stdin");
1375 if (select
[0] == 'y')
1377 pdb_ch
[nch
-1].chainstart
[pdb_ch
[nch
-1].nterpairs
] =
1378 pdba_all
.atom
[i
].resind
- prev_chainstart
;
1379 pdb_ch
[nch
-1].nterpairs
++;
1380 srenew(pdb_ch
[nch
-1].chainstart
,pdb_ch
[nch
-1].nterpairs
+1);
1384 /* set natom for previous chain */
1387 pdb_ch
[nch
-1].natom
=i
-pdb_ch
[nch
-1].start
;
1394 /* check if chain identifier was used before */
1395 for (j
=0; (j
<nch
); j
++)
1397 if (pdb_ch
[j
].chainid
!= ' ' && pdb_ch
[j
].chainid
== ri
->chainid
)
1399 printf("WARNING: Chain identifier '%c' is used in two non-sequential blocks.\n"
1400 "They will be treated as separate chains unless you reorder your file.\n",
1407 srenew(pdb_ch
,maxch
);
1409 pdb_ch
[nch
].chainid
= ri
->chainid
;
1410 pdb_ch
[nch
].chainnum
= ri
->chainnum
;
1411 pdb_ch
[nch
].start
=i
;
1412 pdb_ch
[nch
].bAllWat
=bWat
;
1414 pdb_ch
[nch
].nterpairs
=0;
1416 pdb_ch
[nch
].nterpairs
=1;
1417 snew(pdb_ch
[nch
].chainstart
,pdb_ch
[nch
].nterpairs
+1);
1418 /* modified [nch] to [0] below */
1419 pdb_ch
[nch
].chainstart
[0]=0;
1425 pdb_ch
[nch
-1].natom
=natom
-pdb_ch
[nch
-1].start
;
1427 /* set all the water blocks at the end of the chain */
1428 snew(swap_index
,nch
);
1430 for(i
=0; i
<nch
; i
++)
1431 if (!pdb_ch
[i
].bAllWat
) {
1435 for(i
=0; i
<nch
; i
++)
1436 if (pdb_ch
[i
].bAllWat
) {
1441 printf("Moved all the water blocks to the end\n");
1444 /* copy pdb data and x for all chains */
1445 for (i
=0; (i
<nch
); i
++) {
1447 chains
[i
].chainid
= pdb_ch
[si
].chainid
;
1448 chains
[i
].chainnum
= pdb_ch
[si
].chainnum
;
1449 chains
[i
].bAllWat
= pdb_ch
[si
].bAllWat
;
1450 chains
[i
].nterpairs
= pdb_ch
[si
].nterpairs
;
1451 chains
[i
].chainstart
= pdb_ch
[si
].chainstart
;
1452 snew(chains
[i
].ntdb
,pdb_ch
[si
].nterpairs
);
1453 snew(chains
[i
].ctdb
,pdb_ch
[si
].nterpairs
);
1454 snew(chains
[i
].r_start
,pdb_ch
[si
].nterpairs
);
1455 snew(chains
[i
].r_end
,pdb_ch
[si
].nterpairs
);
1457 snew(chains
[i
].pdba
,1);
1458 init_t_atoms(chains
[i
].pdba
,pdb_ch
[si
].natom
,TRUE
);
1459 snew(chains
[i
].x
,chains
[i
].pdba
->nr
);
1460 for (j
=0; j
<chains
[i
].pdba
->nr
; j
++) {
1461 chains
[i
].pdba
->atom
[j
] = pdba_all
.atom
[pdb_ch
[si
].start
+j
];
1462 snew(chains
[i
].pdba
->atomname
[j
],1);
1463 *chains
[i
].pdba
->atomname
[j
] =
1464 strdup(*pdba_all
.atomname
[pdb_ch
[si
].start
+j
]);
1465 chains
[i
].pdba
->pdbinfo
[j
] = pdba_all
.pdbinfo
[pdb_ch
[si
].start
+j
];
1466 copy_rvec(pdbx
[pdb_ch
[si
].start
+j
],chains
[i
].x
[j
]);
1468 /* Re-index the residues assuming that the indices are continuous */
1469 k
= chains
[i
].pdba
->atom
[0].resind
;
1470 nres
= chains
[i
].pdba
->atom
[chains
[i
].pdba
->nr
-1].resind
- k
+ 1;
1471 chains
[i
].pdba
->nres
= nres
;
1472 for(j
=0; j
< chains
[i
].pdba
->nr
; j
++) {
1473 chains
[i
].pdba
->atom
[j
].resind
-= k
;
1475 srenew(chains
[i
].pdba
->resinfo
,nres
);
1476 for(j
=0; j
<nres
; j
++) {
1477 chains
[i
].pdba
->resinfo
[j
] = pdba_all
.resinfo
[k
+j
];
1478 snew(chains
[i
].pdba
->resinfo
[j
].name
,1);
1479 *chains
[i
].pdba
->resinfo
[j
].name
= strdup(*pdba_all
.resinfo
[k
+j
].name
);
1480 /* make all chain identifiers equal to that of the chain */
1481 chains
[i
].pdba
->resinfo
[j
].chainid
= pdb_ch
[si
].chainid
;
1486 printf("\nMerged %d chains into one molecule definition\n\n",
1487 pdb_ch
[0].nterpairs
);
1489 printf("There are %d chains and %d blocks of water and "
1490 "%d residues with %d atoms\n",
1491 nch
-nwaterchain
,nwaterchain
,
1492 pdba_all
.resinfo
[pdba_all
.atom
[natom
-1].resind
].nr
,natom
);
1494 printf("\n %5s %4s %6s\n","chain","#res","#atoms");
1495 for (i
=0; (i
<nch
); i
++)
1496 printf(" %d '%c' %5d %6d %s\n",
1497 i
+1, chains
[i
].chainid
? chains
[i
].chainid
:'-',
1498 chains
[i
].pdba
->nres
, chains
[i
].pdba
->nr
,
1499 chains
[i
].bAllWat
? "(only water)":"");
1502 check_occupancy(&pdba_all
,opt2fn("-f",NFILE
,fnm
),bVerbose
);
1504 /* Read atomtypes... */
1505 atype
= read_atype(ffdir
,&symtab
);
1507 /* read residue database */
1508 printf("Reading residue database... (%s)\n",forcefield
);
1509 nrtpf
= fflib_search_file_end(ffdir
,".rtp",TRUE
,&rtpf
);
1512 for(i
=0; i
<nrtpf
; i
++) {
1513 read_resall(rtpf
[i
],&nrtp
,&restp
,atype
,&symtab
,FALSE
);
1518 /* Not correct with multiple rtp input files with different bonded types */
1519 fp
=gmx_fio_fopen("new.rtp","w");
1520 print_resall(fp
,nrtp
,restp
,atype
);
1524 /* read hydrogen database */
1525 nah
= read_h_db(ffdir
,&ah
);
1527 /* Read Termini database... */
1528 nNtdb
=read_ter_db(ffdir
,'n',&ntdb
,atype
);
1529 nCtdb
=read_ter_db(ffdir
,'c',&ctdb
,atype
);
1531 top_fn
=ftp2fn(efTOP
,NFILE
,fnm
);
1532 top_file
=gmx_fio_fopen(top_fn
,"w");
1534 sprintf(generator
,"%s - %s",ShortProgram(), GromacsVersion() );
1536 print_top_header(top_file
,top_fn
,generator
,FALSE
,ffdir
,mHmult
);
1543 for(chain
=0; (chain
<nch
); chain
++) {
1544 cc
= &(chains
[chain
]);
1546 /* set pdba, natom and nres to the current chain */
1550 nres
=cc
->pdba
->nres
;
1552 if (cc
->chainid
&& ( cc
->chainid
!= ' ' ) )
1553 printf("Processing chain %d '%c' (%d atoms, %d residues)\n",
1554 chain
+1,cc
->chainid
,natom
,nres
);
1556 printf("Processing chain %d (%d atoms, %d residues)\n",
1557 chain
+1,natom
,nres
);
1559 process_chain(pdba
,x
,bUnA
,bUnA
,bUnA
,bLysMan
,bAspMan
,bGluMan
,
1560 bHisMan
,bArgMan
,bGlnMan
,angle
,distance
,&symtab
,
1561 nrtprename
,rtprename
);
1563 for(i
=0; i
<cc
->nterpairs
; i
++) {
1565 cc
->chainstart
[cc
->nterpairs
] = pdba
->nres
;
1567 find_nc_ter(pdba
,cc
->chainstart
[i
],cc
->chainstart
[i
+1],
1568 &(cc
->r_start
[i
]),&(cc
->r_end
[i
]),rt
);
1571 if ( (cc
->r_start
[i
]<0) || (cc
->r_end
[i
]<0) ) {
1572 printf("Problem with chain definition, or missing terminal residues.\n"
1573 "This chain does not appear to contain a recognized chain molecule.\n"
1574 "If this is incorrect, you can edit residuetypes.dat to modify the behavior.\n");
1581 /* Check for disulfides and other special bonds */
1582 nssbonds
= mk_specbonds(pdba
,x
,bCysMan
,&ssbonds
,bVerbose
);
1584 if (nrtprename
> 0) {
1585 rename_resrtp(pdba
,cc
->nterpairs
,cc
->r_start
,cc
->r_end
,nrtprename
,rtprename
,
1591 sprintf(fn
,"chain.pdb");
1593 sprintf(fn
,"chain_%c%d.pdb",cc
->chainid
,cc
->chainnum
);
1595 write_sto_conf(fn
,title
,pdba
,x
,NULL
,ePBC
,box
);
1599 for(i
=0; i
<cc
->nterpairs
; i
++)
1603 * We first apply a filter so we only have the
1604 * termini that can be applied to the residue in question
1605 * (or a generic terminus if no-residue specific is available).
1607 /* First the N terminus */
1610 tdblist
= filter_ter(nrtp
,restp
,nNtdb
,ntdb
,
1611 *pdba
->resinfo
[cc
->r_start
[i
]].name
,
1612 *pdba
->resinfo
[cc
->r_start
[i
]].rtp
,
1616 printf("No suitable end (N or 5') terminus found in database - assuming this residue\n"
1617 "is already in a terminus-specific form and skipping terminus selection.\n");
1622 if(bTerMan
&& ntdblist
>1)
1624 cc
->ntdb
[i
] = choose_ter(ntdblist
,tdblist
,"Select start terminus type");
1628 cc
->ntdb
[i
] = tdblist
[0];
1631 printf("Start terminus: %s\n",(cc
->ntdb
[i
])->name
);
1640 /* And the C terminus */
1643 tdblist
= filter_ter(nrtp
,restp
,nCtdb
,ctdb
,
1644 *pdba
->resinfo
[cc
->r_end
[i
]].name
,
1645 *pdba
->resinfo
[cc
->r_end
[i
]].rtp
,
1649 printf("No suitable end (C or 3') terminus found in database - assuming this residue\n"
1650 "is already in a terminus-specific form and skipping terminus selection.\n");
1655 if(bTerMan
&& ntdblist
>1)
1657 cc
->ctdb
[i
] = choose_ter(ntdblist
,tdblist
,"Select end terminus type");
1661 cc
->ctdb
[i
] = tdblist
[0];
1663 printf("End terminus: %s\n",(cc
->ctdb
[i
])->name
);
1672 /* lookup hackblocks and rtp for all residues */
1673 get_hackblocks_rtp(&hb_chain
, &restp_chain
,
1674 nrtp
, restp
, pdba
->nres
, pdba
->resinfo
,
1675 cc
->nterpairs
, cc
->ntdb
, cc
->ctdb
, cc
->r_start
, cc
->r_end
);
1676 /* ideally, now we would not need the rtp itself anymore, but do
1677 everything using the hb and restp arrays. Unfortunately, that
1678 requires some re-thinking of code in gen_vsite.c, which I won't
1679 do now :( AF 26-7-99 */
1681 rename_atoms(NULL
,ffdir
,
1682 pdba
,&symtab
,restp_chain
,FALSE
,rt
,FALSE
,bVerbose
);
1684 match_atomnames_with_rtp(restp_chain
,hb_chain
,pdba
,x
,bVerbose
);
1687 block
= new_blocka();
1689 sort_pdbatoms(pdba
->nres
,restp_chain
,hb_chain
,
1690 natom
,&pdba
,&x
,block
,&gnames
);
1691 natom
= remove_duplicate_atoms(pdba
,x
,bVerbose
);
1692 if (ftp2bSet(efNDX
,NFILE
,fnm
)) {
1694 fprintf(stderr
,"WARNING: with the -remh option the generated "
1695 "index file (%s) might be useless\n"
1696 "(the index file is generated before hydrogens are added)",
1697 ftp2fn(efNDX
,NFILE
,fnm
));
1699 write_index(ftp2fn(efNDX
,NFILE
,fnm
),block
,gnames
);
1701 for(i
=0; i
< block
->nr
; i
++)
1706 fprintf(stderr
,"WARNING: "
1707 "without sorting no check for duplicate atoms can be done\n");
1710 /* Generate Hydrogen atoms (and termini) in the sequence */
1711 natom
=add_h(&pdba
,&x
,nah
,ah
,
1712 cc
->nterpairs
,cc
->ntdb
,cc
->ctdb
,cc
->r_start
,cc
->r_end
,bAllowMissing
,
1713 NULL
,NULL
,TRUE
,FALSE
);
1714 printf("Now there are %d residues with %d atoms\n",
1715 pdba
->nres
,pdba
->nr
);
1716 if (debug
) write_pdbfile(debug
,title
,pdba
,x
,ePBC
,box
,' ',0,NULL
,TRUE
);
1719 for(i
=0; (i
<natom
); i
++)
1720 fprintf(debug
,"Res %s%d atom %d %s\n",
1721 *(pdba
->resinfo
[pdba
->atom
[i
].resind
].name
),
1722 pdba
->resinfo
[pdba
->atom
[i
].resind
].nr
,i
+1,*pdba
->atomname
[i
]);
1724 strcpy(posre_fn
,ftp2fn(efITP
,NFILE
,fnm
));
1726 /* make up molecule name(s) */
1728 k
= (cc
->nterpairs
>0 && cc
->r_start
[0]>=0) ? cc
->r_start
[0] : 0;
1730 gmx_residuetype_get_type(rt
,*pdba
->resinfo
[k
].name
,&p_restype
);
1736 sprintf(molname
,"Water");
1740 this_chainid
= cc
->chainid
;
1742 /* Add the chain id if we have one */
1743 if(this_chainid
!= ' ')
1745 sprintf(buf
,"_chain_%c",this_chainid
);
1749 /* Check if there have been previous chains with the same id */
1751 for(k
=0;k
<chain
;k
++)
1753 if(cc
->chainid
== chains
[k
].chainid
)
1758 /* Add the number for this chain identifier if there are multiple copies */
1762 sprintf(buf
,"%d",nid_used
+1);
1766 if(strlen(suffix
)>0)
1768 sprintf(molname
,"%s%s",p_restype
,suffix
);
1772 strcpy(molname
,p_restype
);
1776 if ((nch
-nwaterchain
>1) && !cc
->bAllWat
) {
1778 strcpy(itp_fn
,top_fn
);
1779 printf("Chain time...\n");
1780 c
=strrchr(itp_fn
,'.');
1781 sprintf(c
,"_%s.itp",molname
);
1782 c
=strrchr(posre_fn
,'.');
1783 sprintf(c
,"_%s.itp",molname
);
1784 if (strcmp(itp_fn
,posre_fn
) == 0) {
1785 strcpy(buf_fn
,posre_fn
);
1786 c
= strrchr(buf_fn
,'.');
1788 sprintf(posre_fn
,"%s_pr.itp",buf_fn
);
1792 srenew(incls
,nincl
);
1793 incls
[nincl
-1]=strdup(itp_fn
);
1794 itp_file
=gmx_fio_fopen(itp_fn
,"w");
1798 srenew(mols
,nmol
+1);
1800 mols
[nmol
].name
= strdup("SOL");
1801 mols
[nmol
].nr
= pdba
->nres
;
1803 mols
[nmol
].name
= strdup(molname
);
1809 print_top_comment(itp_file
,itp_fn
,generator
,ffdir
,TRUE
);
1819 pdb2top(top_file2
,posre_fn
,molname
,pdba
,&x
,atype
,&symtab
,
1821 restp_chain
,hb_chain
,
1822 cc
->nterpairs
,cc
->ntdb
,cc
->ctdb
,cc
->r_start
,cc
->r_end
,bAllowMissing
,
1823 bVsites
,bVsiteAromatics
,forcefield
,ffdir
,
1824 mHmult
,nssbonds
,ssbonds
,
1825 long_bond_dist
,short_bond_dist
,bDeuterate
,bChargeGroups
,bCmap
,
1826 bRenumRes
,bRTPresname
);
1829 write_posres(posre_fn
,pdba
,posre_fc
);
1832 gmx_fio_fclose(itp_file
);
1834 /* pdba and natom have been reassigned somewhere so: */
1839 if (cc
->chainid
== ' ')
1840 sprintf(fn
,"chain.pdb");
1842 sprintf(fn
,"chain_%c.pdb",cc
->chainid
);
1843 cool_quote(quote
,255,NULL
);
1844 write_sto_conf(fn
,quote
,pdba
,x
,NULL
,ePBC
,box
);
1848 if (watermodel
== NULL
) {
1849 for(chain
=0; chain
<nch
; chain
++) {
1850 if (chains
[chain
].bAllWat
) {
1851 gmx_fatal(FARGS
,"You have chosen not to include a water model, but there is water in the input file. Select a water model or remove the water from your input file.");
1855 sprintf(buf_fn
,"%s%c%s.itp",ffdir
,DIR_SEPARATOR
,watermodel
);
1856 if (!fflib_fexist(buf_fn
)) {
1857 gmx_fatal(FARGS
,"The topology file '%s' for the selected water model '%s' can not be found in the force field directory. Select a different water model.",
1862 print_top_mols(top_file
,title
,ffdir
,watermodel
,nincl
,incls
,nmol
,mols
);
1863 gmx_fio_fclose(top_file
);
1865 gmx_residuetype_destroy(rt
);
1867 /* now merge all chains back together */
1870 for (i
=0; (i
<nch
); i
++) {
1871 natom
+=chains
[i
].pdba
->nr
;
1872 nres
+=chains
[i
].pdba
->nres
;
1875 init_t_atoms(atoms
,natom
,FALSE
);
1876 for(i
=0; i
< atoms
->nres
; i
++)
1877 sfree(atoms
->resinfo
[i
].name
);
1878 sfree(atoms
->resinfo
);
1880 snew(atoms
->resinfo
,nres
);
1884 for (i
=0; (i
<nch
); i
++) {
1886 printf("Including chain %d in system: %d atoms %d residues\n",
1887 i
+1,chains
[i
].pdba
->nr
,chains
[i
].pdba
->nres
);
1888 for (j
=0; (j
<chains
[i
].pdba
->nr
); j
++) {
1889 atoms
->atom
[k
]=chains
[i
].pdba
->atom
[j
];
1890 atoms
->atom
[k
].resind
+= l
; /* l is processed nr of residues */
1891 atoms
->atomname
[k
]=chains
[i
].pdba
->atomname
[j
];
1892 atoms
->resinfo
[atoms
->atom
[k
].resind
].chainid
= chains
[i
].chainid
;
1893 copy_rvec(chains
[i
].x
[j
],x
[k
]);
1896 for (j
=0; (j
<chains
[i
].pdba
->nres
); j
++) {
1897 atoms
->resinfo
[l
] = chains
[i
].pdba
->resinfo
[j
];
1899 atoms
->resinfo
[l
].name
= atoms
->resinfo
[l
].rtp
;
1906 fprintf(stderr
,"Now there are %d atoms and %d residues\n",k
,l
);
1907 print_sums(atoms
, TRUE
);
1910 fprintf(stderr
,"\nWriting coordinate file...\n");
1911 clear_rvec(box_space
);
1913 gen_box(0,atoms
->nr
,x
,box
,box_space
,FALSE
);
1914 write_sto_conf(ftp2fn(efSTO
,NFILE
,fnm
),title
,atoms
,x
,NULL
,ePBC
,box
);
1916 printf("\t\t--------- PLEASE NOTE ------------\n");
1917 printf("You have successfully generated a topology from: %s.\n",
1918 opt2fn("-f",NFILE
,fnm
));
1919 if (watermodel
!= NULL
) {
1920 printf("The %s force field and the %s water model are used.\n",
1923 printf("The %s force field is used.\n",
1926 printf("\t\t--------- ETON ESAELP ------------\n");