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.
48 #include "gromacs/fileio/confio.h"
51 #include "gromacs/fileio/futil.h"
52 #include "gmx_fatal.h"
63 static void copy_atom(t_atoms *atoms1, int a1, t_atoms *atoms2, int a2)
65 atoms2->atom[a2] = atoms1->atom[a1];
66 snew(atoms2->atomname[a2], 1);
67 *atoms2->atomname[a2] = strdup(*atoms1->atomname[a1]);
70 static atom_id pdbasearch_atom(const char *name, int resind, t_atoms *pdba,
71 const char *searchtype, gmx_bool bAllowMissing)
75 for (i = 0; (i < pdba->nr) && (pdba->atom[i].resind != resind); i++)
80 return search_atom(name, i, pdba,
81 searchtype, bAllowMissing);
84 static void hacksearch_atom(int *ii, int *jj, char *name,
85 int nab[], t_hack *ab[],
86 int resind, t_atoms *pdba)
96 for (i = 0; (i < pdba->nr) && (pdba->atom[i].resind != resind); i++)
100 for (; (i < pdba->nr) && (pdba->atom[i].resind == resind) && (*ii < 0); i++)
102 for (j = 0; (j < nab[i]) && (*ii < 0); j++)
104 if (ab[i][j].nname && strcmp(name, ab[i][j].nname) == 0)
115 void dump_ab(FILE *out, int natom, int nab[], t_hack *ab[], gmx_bool bHeader)
119 #define SS(s) (s) ? (s) : "-"
123 fprintf(out, "ADDBLOCK (t_hack) natom=%d\n"
124 "%4s %2s %-4s %-4s %2s %-4s %-4s %-4s %-4s %1s %s\n",
125 natom, "atom", "nr", "old", "new", "tp", "ai", "aj", "ak", "al", "a", "x");
127 for (i = 0; i < natom; i++)
129 for (j = 0; j < nab[i]; j++)
131 fprintf(out, "%4d %2d %-4s %-4s %2d %-4s %-4s %-4s %-4s %s %g %g %g\n",
132 i+1, ab[i][j].nr, SS(ab[i][j].oname), SS(ab[i][j].nname),
134 SS(ab[i][j].AI), SS(ab[i][j].AJ),
135 SS(ab[i][j].AK), SS(ab[i][j].AL),
136 ab[i][j].atom ? "+" : "",
137 ab[i][j].newx[XX], ab[i][j].newx[YY], ab[i][j].newx[ZZ]);
143 static t_hackblock *get_hackblocks(t_atoms *pdba, int nah, t_hackblock ah[],
145 t_hackblock **ntdb, t_hackblock **ctdb,
149 t_hackblock *hb, *ahptr;
152 snew(hb, pdba->nres);
153 /* first the termini */
154 for (i = 0; i < nterpairs; i++)
158 copy_t_hackblock(ntdb[i], &hb[rN[i]]);
162 merge_t_hackblock(ctdb[i], &hb[rC[i]]);
165 /* then the whole hdb */
166 for (rnr = 0; rnr < pdba->nres; rnr++)
168 ahptr = search_h_db(nah, ah, *pdba->resinfo[rnr].rtp);
171 if (hb[rnr].name == NULL)
173 hb[rnr].name = strdup(ahptr->name);
175 merge_hacks(ahptr, &hb[rnr]);
181 static void expand_hackblocks_one(t_hackblock *hbr, char *atomname,
182 int *nabi, t_hack **abi, gmx_bool bN, gmx_bool bC)
187 /* we'll recursively add atoms to atoms */
188 for (j = 0; j < hbr->nhack; j++)
190 /* first check if we're in the N- or C-terminus, then we should ignore
191 all hacks involving atoms from resp. previous or next residue
192 (i.e. which name begins with '-' (N) or '+' (C) */
194 if (bN) /* N-terminus: ignore '-' */
196 for (k = 0; k < 4 && hbr->hack[j].a[k] && !bIgnore; k++)
198 bIgnore = hbr->hack[j].a[k][0] == '-';
201 if (bC) /* C-terminus: ignore '+' */
203 for (k = 0; k < 4 && hbr->hack[j].a[k] && !bIgnore; k++)
205 bIgnore = hbr->hack[j].a[k][0] == '+';
208 /* must be either hdb entry (tp>0) or add from tdb (oname==NULL)
209 and first control aton (AI) matches this atom or
210 delete/replace from tdb (oname!=NULL) and oname matches this atom */
213 fprintf(debug, " %s", hbr->hack[j].oname ? hbr->hack[j].oname : hbr->hack[j].AI);
217 ( ( ( hbr->hack[j].tp > 0 || hbr->hack[j].oname == NULL ) &&
218 strcmp(atomname, hbr->hack[j].AI) == 0 ) ||
219 ( hbr->hack[j].oname != NULL &&
220 strcmp(atomname, hbr->hack[j].oname) == 0) ) )
222 /* now expand all hacks for this atom */
225 fprintf(debug, " +%dh", hbr->hack[j].nr);
227 srenew(*abi, *nabi + hbr->hack[j].nr);
228 for (k = 0; k < hbr->hack[j].nr; k++)
230 copy_t_hack(&hbr->hack[j], &(*abi)[*nabi + k]);
231 (*abi)[*nabi + k].bXSet = FALSE;
232 /* if we're adding (oname==NULL) and don't have a new name (nname)
233 yet, build it from atomname */
234 if ( (*abi)[*nabi + k].nname == NULL)
236 if ( (*abi)[*nabi + k].oname == NULL)
238 (*abi)[*nabi + k].nname = strdup(atomname);
239 (*abi)[*nabi + k].nname[0] = 'H';
246 fprintf(debug, "Hack '%s' %d, replacing nname '%s' with '%s' (old name '%s')\n",
248 (*abi)[*nabi + k].nname, hbr->hack[j].nname,
249 (*abi)[*nabi + k].oname ? (*abi)[*nabi + k].oname : "");
251 sfree((*abi)[*nabi + k].nname);
252 (*abi)[*nabi + k].nname = strdup(hbr->hack[j].nname);
255 if (hbr->hack[j].tp == 10 && k == 2)
257 /* This is a water virtual site, not a hydrogen */
258 /* Ugly hardcoded name hack */
259 (*abi)[*nabi + k].nname[0] = 'M';
261 else if (hbr->hack[j].tp == 11 && k >= 2)
263 /* This is a water lone pair, not a hydrogen */
264 /* Ugly hardcoded name hack */
265 srenew((*abi)[*nabi + k].nname, 4);
266 (*abi)[*nabi + k].nname[0] = 'L';
267 (*abi)[*nabi + k].nname[1] = 'P';
268 (*abi)[*nabi + k].nname[2] = '1' + k - 2;
269 (*abi)[*nabi + k].nname[3] = '\0';
271 else if (hbr->hack[j].nr > 1)
273 /* adding more than one atom, number them */
274 l = strlen((*abi)[*nabi + k].nname);
275 srenew((*abi)[*nabi + k].nname, l+2);
276 (*abi)[*nabi + k].nname[l] = '1' + k;
277 (*abi)[*nabi + k].nname[l+1] = '\0';
280 (*nabi) += hbr->hack[j].nr;
282 /* add hacks to atoms we've just added */
283 if (hbr->hack[j].tp > 0 || hbr->hack[j].oname == NULL)
285 for (k = 0; k < hbr->hack[j].nr; k++)
287 expand_hackblocks_one(hbr, (*abi)[*nabi-hbr->hack[j].nr+k].nname,
295 static void expand_hackblocks(t_atoms *pdba, t_hackblock hb[],
296 int nab[], t_hack *ab[],
297 int nterpairs, int *rN, int *rC)
302 for (i = 0; i < pdba->nr; i++)
305 for (j = 0; j < nterpairs && !bN; j++)
307 bN = pdba->atom[i].resind == rN[j];
310 for (j = 0; j < nterpairs && !bC; j++)
312 bC = pdba->atom[i].resind == rC[j];
315 /* add hacks to this atom */
316 expand_hackblocks_one(&hb[pdba->atom[i].resind], *pdba->atomname[i],
317 &nab[i], &ab[i], bN, bC);
321 fprintf(debug, "\n");
325 static int check_atoms_present(t_atoms *pdba, int nab[], t_hack *ab[])
327 int i, j, k, d, rnr, nadd;
330 for (i = 0; i < pdba->nr; i++)
332 rnr = pdba->atom[i].resind;
333 for (j = 0; j < nab[i]; j++)
335 if (ab[i][j].oname == NULL)
338 if (ab[i][j].nname == NULL)
340 gmx_incons("ab[i][j].nname not allocated");
342 /* check if the atom is already present */
343 k = pdbasearch_atom(ab[i][j].nname, rnr, pdba, "check", TRUE);
346 /* We found the added atom. */
347 ab[i][j].bAlreadyPresent = TRUE;
350 fprintf(debug, "Atom '%s' in residue '%s' %d is already present\n",
352 *pdba->resinfo[rnr].name, pdba->resinfo[rnr].nr);
357 ab[i][j].bAlreadyPresent = FALSE;
358 /* count how many atoms we'll add */
362 else if (ab[i][j].nname == NULL)
373 static void calc_all_pos(t_atoms *pdba, rvec x[], int nab[], t_hack *ab[],
374 gmx_bool bCheckMissing)
376 int i, j, ii, jj, m, ia, d, rnr, l = 0;
378 rvec xa[4]; /* control atoms for calc_h_pos */
379 rvec xh[MAXH]; /* hydrogen positions from calc_h_pos */
384 for (i = 0; i < pdba->nr; i++)
386 rnr = pdba->atom[i].resind;
387 for (j = 0; j < nab[i]; j += ab[i][j].nr)
389 /* check if we're adding: */
390 if (ab[i][j].oname == NULL && ab[i][j].tp > 0)
393 for (m = 0; (m < ab[i][j].nctl && bFoundAll); m++)
395 ia = pdbasearch_atom(ab[i][j].a[m], rnr, pdba,
396 bCheckMissing ? "atom" : "check",
400 /* not found in original atoms, might still be in t_hack (ab) */
401 hacksearch_atom(&ii, &jj, ab[i][j].a[m], nab, ab, rnr, pdba);
404 copy_rvec(ab[ii][jj].newx, xa[m]);
411 gmx_fatal(FARGS, "Atom %s not found in residue %s %d"
413 " while adding hydrogens",
415 *pdba->resinfo[rnr].name,
416 pdba->resinfo[rnr].nr,
417 *pdba->resinfo[rnr].rtp);
423 copy_rvec(x[ia], xa[m]);
428 for (m = 0; (m < MAXH); m++)
430 for (d = 0; d < DIM; d++)
442 calc_h_pos(ab[i][j].tp, xa, xh, &l);
443 for (m = 0; m < ab[i][j].nr; m++)
445 copy_rvec(xh[m], ab[i][j+m].newx);
446 ab[i][j+m].bXSet = TRUE;
454 static void free_ab(int natoms, int *nab, t_hack **ab)
458 for (i = 0; i < natoms; i++)
460 free_t_hack(nab[i], &ab[i]);
466 static int add_h_low(t_atoms **pdbaptr, rvec *xptr[],
467 int nah, t_hackblock ah[],
468 int nterpairs, t_hackblock **ntdb, t_hackblock **ctdb,
469 int *rN, int *rC, gmx_bool bCheckMissing,
470 int **nabptr, t_hack ***abptr,
471 gmx_bool bUpdate_pdba, gmx_bool bKeep_old_pdba)
473 t_atoms *newpdba = NULL, *pdba = NULL;
475 int i, newi, j, d, natoms, nalreadypresent;
482 /* set flags for adding hydrogens (according to hdb) */
488 /* the first time these will be pointers to NULL, but we will
489 return in them the completed arrays, which we will get back
496 fprintf(debug, "pointer to ab found\n");
506 /* WOW, everything was already figured out */
507 bUpdate_pdba = FALSE;
510 fprintf(debug, "pointer to non-null ab found\n");
515 /* We'll have to do all the hard work */
517 /* first get all the hackblocks for each residue: */
518 hb = get_hackblocks(pdba, nah, ah, nterpairs, ntdb, ctdb, rN, rC);
521 dump_hb(debug, pdba->nres, hb);
524 /* expand the hackblocks to atom level */
527 expand_hackblocks(pdba, hb, nab, ab, nterpairs, rN, rC);
528 free_t_hackblock(pdba->nres, &hb);
533 fprintf(debug, "before calc_all_pos\n");
534 dump_ab(debug, natoms, nab, ab, TRUE);
537 /* Now calc the positions */
538 calc_all_pos(pdba, *xptr, nab, ab, bCheckMissing);
542 fprintf(debug, "after calc_all_pos\n");
543 dump_ab(debug, natoms, nab, ab, TRUE);
548 /* we don't have to add atoms that are already present in pdba,
549 so we will remove them from the ab (t_hack) */
550 nadd = check_atoms_present(pdba, nab, ab);
553 fprintf(debug, "removed add hacks that were already in pdba:\n");
554 dump_ab(debug, natoms, nab, ab, TRUE);
555 fprintf(debug, "will be adding %d atoms\n", nadd);
558 /* Copy old atoms, making space for new ones */
560 init_t_atoms(newpdba, natoms+nadd, FALSE);
561 newpdba->nres = pdba->nres;
562 sfree(newpdba->resinfo);
563 newpdba->resinfo = pdba->resinfo;
571 fprintf(debug, "snew xn for %d old + %d new atoms %d total)\n",
572 natoms, nadd, natoms+nadd);
577 /* There is nothing to do: return now */
580 free_ab(natoms, nab, ab);
586 snew(xn, natoms+nadd);
588 for (i = 0; (i < natoms); i++)
590 /* check if this atom wasn't scheduled for deletion */
591 if (nab[i] == 0 || (ab[i][0].nname != NULL) )
593 if (newi >= natoms+nadd)
595 /*gmx_fatal(FARGS,"Not enough space for adding atoms");*/
597 srenew(xn, natoms+nadd);
600 srenew(newpdba->atom, natoms+nadd);
601 srenew(newpdba->atomname, natoms+nadd);
607 fprintf(debug, "(%3d) %3d %4s %4s%3d %3d",
608 i+1, newi+1, *pdba->atomname[i],
609 *pdba->resinfo[pdba->atom[i].resind].name,
610 pdba->resinfo[pdba->atom[i].resind].nr, nab[i]);
614 copy_atom(pdba, i, newpdba, newi);
616 copy_rvec((*xptr)[i], xn[newi]);
617 /* process the hacks for this atom */
619 for (j = 0; j < nab[i]; j++)
621 if (ab[i][j].oname == NULL) /* add */
624 if (newi >= natoms+nadd)
626 /* gmx_fatal(FARGS,"Not enough space for adding atoms");*/
628 srenew(xn, natoms+nadd);
631 srenew(newpdba->atom, natoms+nadd);
632 srenew(newpdba->atomname, natoms+nadd);
638 newpdba->atom[newi].resind = pdba->atom[i].resind;
642 fprintf(debug, " + %d", newi+1);
645 if (ab[i][j].nname != NULL &&
646 (ab[i][j].oname == NULL ||
647 strcmp(ab[i][j].oname, *newpdba->atomname[newi]) == 0))
650 if (ab[i][j].oname == NULL && ab[i][j].bAlreadyPresent)
652 /* This atom is already present, copy it from the input. */
656 copy_atom(pdba, i+nalreadypresent, newpdba, newi);
658 copy_rvec((*xptr)[i+nalreadypresent], xn[newi]);
666 fprintf(debug, "Replacing %d '%s' with (old name '%s') %s\n",
668 (newpdba->atomname[newi] && *newpdba->atomname[newi]) ? *newpdba->atomname[newi] : "",
669 ab[i][j].oname ? ab[i][j].oname : "",
672 snew(newpdba->atomname[newi], 1);
673 *newpdba->atomname[newi] = strdup(ab[i][j].nname);
674 if (ab[i][j].oname != NULL && ab[i][j].atom) /* replace */
675 { /* newpdba->atom[newi].m = ab[i][j].atom->m; */
676 /* newpdba->atom[newi].q = ab[i][j].atom->q; */
677 /* newpdba->atom[newi].type = ab[i][j].atom->type; */
682 copy_rvec(ab[i][j].newx, xn[newi]);
685 if (bUpdate_pdba && debug)
687 fprintf(debug, " %s %g %g", *newpdba->atomname[newi],
688 newpdba->atom[newi].m, newpdba->atom[newi].q);
693 i += nalreadypresent;
696 fprintf(debug, "\n");
713 free_ab(natoms, nab, ab);
720 for (i = 0; i < natoms; i++)
722 /* Do not free the atomname string itself, it might be in symtab */
723 /* sfree(*(pdba->atomname[i])); */
724 /* sfree(pdba->atomname[i]); */
726 sfree(pdba->atomname);
728 sfree(pdba->pdbinfo);
744 void deprotonate(t_atoms *atoms, rvec *x)
749 for (i = 0; i < atoms->nr; i++)
751 if ( (*atoms->atomname[i])[0] != 'H')
753 atoms->atomname[j] = atoms->atomname[i];
754 atoms->atom[j] = atoms->atom[i];
755 copy_rvec(x[i], x[j]);
762 int add_h(t_atoms **pdbaptr, rvec *xptr[],
763 int nah, t_hackblock ah[],
764 int nterpairs, t_hackblock **ntdb, t_hackblock **ctdb,
765 int *rN, int *rC, gmx_bool bAllowMissing,
766 int **nabptr, t_hack ***abptr,
767 gmx_bool bUpdate_pdba, gmx_bool bKeep_old_pdba)
769 int nold, nnew, niter;
771 /* Here we loop to be able to add atoms to added atoms.
772 * We should not check for missing atoms here.
779 nnew = add_h_low(pdbaptr, xptr, nah, ah, nterpairs, ntdb, ctdb, rN, rC, FALSE,
780 nabptr, abptr, bUpdate_pdba, bKeep_old_pdba);
784 gmx_fatal(FARGS, "More than 100 iterations of add_h. Maybe you are trying to replace an added atom (this is not supported)?");
791 /* Call add_h_low once more, now only for the missing atoms check */
792 add_h_low(pdbaptr, xptr, nah, ah, nterpairs, ntdb, ctdb, rN, rC, TRUE,
793 nabptr, abptr, bUpdate_pdba, bKeep_old_pdba);
799 int protonate(t_atoms **atomsptr, rvec **xptr, t_protonate *protdata)
803 gmx_bool bUpdate_pdba, bKeep_old_pdba;
804 int nntdb, nctdb, nt, ct;
808 if (!protdata->bInit)
812 fprintf(debug, "protonate: Initializing protdata\n");
815 /* set forcefield to use: */
816 strcpy(protdata->FF, "oplsaa.ff");
818 /* get the databases: */
819 protdata->nah = read_h_db(protdata->FF, &protdata->ah);
820 open_symtab(&protdata->tab);
821 protdata->atype = read_atype(protdata->FF, &protdata->tab);
822 nntdb = read_ter_db(protdata->FF, 'n', &protdata->ntdb, protdata->atype);
825 gmx_fatal(FARGS, "no N-terminus database");
827 nctdb = read_ter_db(protdata->FF, 'c', &protdata->ctdb, protdata->atype);
830 gmx_fatal(FARGS, "no C-terminus database");
833 /* set terminus types: -NH3+ (different for Proline) and -COO- */
835 snew(protdata->sel_ntdb, NTERPAIRS);
836 snew(protdata->sel_ctdb, NTERPAIRS);
838 if (nntdb >= 4 && nctdb >= 2)
840 /* Yuk, yuk, hardcoded default termini selections !!! */
841 if (strncmp(*atoms->resinfo[atoms->atom[atoms->nr-1].resind].name, "PRO", 3) == 0)
856 protdata->sel_ntdb[0] = &(protdata->ntdb[nt]);
857 protdata->sel_ctdb[0] = &(protdata->ctdb[ct]);
859 /* set terminal residue numbers: */
860 snew(protdata->rN, NTERPAIRS);
861 snew(protdata->rC, NTERPAIRS);
863 protdata->rC[0] = atoms->atom[atoms->nr-1].resind;
865 /* keep unprotonated topology: */
866 protdata->upatoms = atoms;
867 /* we don't have these yet: */
868 protdata->patoms = NULL;
870 bKeep_old_pdba = TRUE;
872 /* clear hackblocks: */
873 protdata->nab = NULL;
876 /* set flag to show we're initialized: */
877 protdata->bInit = TRUE;
883 fprintf(debug, "protonate: using available protdata\n");
885 /* add_h will need the unprotonated topology again: */
886 atoms = protdata->upatoms;
887 bUpdate_pdba = FALSE;
888 bKeep_old_pdba = FALSE;
892 nadd = add_h(&atoms, xptr, protdata->nah, protdata->ah,
893 NTERPAIRS, protdata->sel_ntdb, protdata->sel_ctdb,
894 protdata->rN, protdata->rC, TRUE,
895 &protdata->nab, &protdata->ab, bUpdate_pdba, bKeep_old_pdba);
896 if (!protdata->patoms)
898 /* store protonated topology */
899 protdata->patoms = atoms;
901 *atomsptr = protdata->patoms;
904 fprintf(debug, "natoms: %d -> %d (nadd=%d)\n",
905 protdata->upatoms->nr, protdata->patoms->nr, nadd);