2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, 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.
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/gmxpreprocess/calch.h"
46 #include "gromacs/gmxpreprocess/h_db.h"
47 #include "gromacs/gmxpreprocess/pgutil.h"
48 #include "gromacs/gmxpreprocess/resall.h"
49 #include "gromacs/gmxpreprocess/ter_db.h"
50 #include "gromacs/legacyheaders/network.h"
51 #include "gromacs/legacyheaders/typedefs.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/topology/symtab.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/futil.h"
57 #include "gromacs/utility/smalloc.h"
59 static void copy_atom(t_atoms *atoms1, int a1, t_atoms *atoms2, int a2)
61 atoms2->atom[a2] = atoms1->atom[a1];
62 snew(atoms2->atomname[a2], 1);
63 *atoms2->atomname[a2] = gmx_strdup(*atoms1->atomname[a1]);
66 static atom_id pdbasearch_atom(const char *name, int resind, t_atoms *pdba,
67 const char *searchtype, gmx_bool bAllowMissing)
71 for (i = 0; (i < pdba->nr) && (pdba->atom[i].resind != resind); i++)
76 return search_atom(name, i, pdba,
77 searchtype, bAllowMissing);
80 static void hacksearch_atom(int *ii, int *jj, char *name,
81 int nab[], t_hack *ab[],
82 int resind, t_atoms *pdba)
92 for (i = 0; (i < pdba->nr) && (pdba->atom[i].resind != resind); i++)
96 for (; (i < pdba->nr) && (pdba->atom[i].resind == resind) && (*ii < 0); i++)
98 for (j = 0; (j < nab[i]) && (*ii < 0); j++)
100 if (ab[i][j].nname && strcmp(name, ab[i][j].nname) == 0)
111 void dump_ab(FILE *out, int natom, int nab[], t_hack *ab[], gmx_bool bHeader)
115 #define SS(s) (s) ? (s) : "-"
119 fprintf(out, "ADDBLOCK (t_hack) natom=%d\n"
120 "%4s %2s %-4s %-4s %2s %-4s %-4s %-4s %-4s %1s %s\n",
121 natom, "atom", "nr", "old", "new", "tp", "ai", "aj", "ak", "al", "a", "x");
123 for (i = 0; i < natom; i++)
125 for (j = 0; j < nab[i]; j++)
127 fprintf(out, "%4d %2d %-4s %-4s %2d %-4s %-4s %-4s %-4s %s %g %g %g\n",
128 i+1, ab[i][j].nr, SS(ab[i][j].oname), SS(ab[i][j].nname),
130 SS(ab[i][j].AI), SS(ab[i][j].AJ),
131 SS(ab[i][j].AK), SS(ab[i][j].AL),
132 ab[i][j].atom ? "+" : "",
133 ab[i][j].newx[XX], ab[i][j].newx[YY], ab[i][j].newx[ZZ]);
139 static t_hackblock *get_hackblocks(t_atoms *pdba, int nah, t_hackblock ah[],
141 t_hackblock **ntdb, t_hackblock **ctdb,
145 t_hackblock *hb, *ahptr;
148 snew(hb, pdba->nres);
149 /* first the termini */
150 for (i = 0; i < nterpairs; i++)
154 copy_t_hackblock(ntdb[i], &hb[rN[i]]);
158 merge_t_hackblock(ctdb[i], &hb[rC[i]]);
161 /* then the whole hdb */
162 for (rnr = 0; rnr < pdba->nres; rnr++)
164 ahptr = search_h_db(nah, ah, *pdba->resinfo[rnr].rtp);
167 if (hb[rnr].name == NULL)
169 hb[rnr].name = gmx_strdup(ahptr->name);
171 merge_hacks(ahptr, &hb[rnr]);
177 static void expand_hackblocks_one(t_hackblock *hbr, char *atomname,
178 int *nabi, t_hack **abi, gmx_bool bN, gmx_bool bC)
183 /* we'll recursively add atoms to atoms */
184 for (j = 0; j < hbr->nhack; j++)
186 /* first check if we're in the N- or C-terminus, then we should ignore
187 all hacks involving atoms from resp. previous or next residue
188 (i.e. which name begins with '-' (N) or '+' (C) */
190 if (bN) /* N-terminus: ignore '-' */
192 for (k = 0; k < 4 && hbr->hack[j].a[k] && !bIgnore; k++)
194 bIgnore = hbr->hack[j].a[k][0] == '-';
197 if (bC) /* C-terminus: ignore '+' */
199 for (k = 0; k < 4 && hbr->hack[j].a[k] && !bIgnore; k++)
201 bIgnore = hbr->hack[j].a[k][0] == '+';
204 /* must be either hdb entry (tp>0) or add from tdb (oname==NULL)
205 and first control aton (AI) matches this atom or
206 delete/replace from tdb (oname!=NULL) and oname matches this atom */
209 fprintf(debug, " %s", hbr->hack[j].oname ? hbr->hack[j].oname : hbr->hack[j].AI);
213 ( ( ( hbr->hack[j].tp > 0 || hbr->hack[j].oname == NULL ) &&
214 strcmp(atomname, hbr->hack[j].AI) == 0 ) ||
215 ( hbr->hack[j].oname != NULL &&
216 strcmp(atomname, hbr->hack[j].oname) == 0) ) )
218 /* now expand all hacks for this atom */
221 fprintf(debug, " +%dh", hbr->hack[j].nr);
223 srenew(*abi, *nabi + hbr->hack[j].nr);
224 for (k = 0; k < hbr->hack[j].nr; k++)
226 copy_t_hack(&hbr->hack[j], &(*abi)[*nabi + k]);
227 (*abi)[*nabi + k].bXSet = FALSE;
228 /* if we're adding (oname==NULL) and don't have a new name (nname)
229 yet, build it from atomname */
230 if ( (*abi)[*nabi + k].nname == NULL)
232 if ( (*abi)[*nabi + k].oname == NULL)
234 (*abi)[*nabi + k].nname = gmx_strdup(atomname);
235 (*abi)[*nabi + k].nname[0] = 'H';
242 fprintf(debug, "Hack '%s' %d, replacing nname '%s' with '%s' (old name '%s')\n",
244 (*abi)[*nabi + k].nname, hbr->hack[j].nname,
245 (*abi)[*nabi + k].oname ? (*abi)[*nabi + k].oname : "");
247 sfree((*abi)[*nabi + k].nname);
248 (*abi)[*nabi + k].nname = gmx_strdup(hbr->hack[j].nname);
251 if (hbr->hack[j].tp == 10 && k == 2)
253 /* This is a water virtual site, not a hydrogen */
254 /* Ugly hardcoded name hack */
255 (*abi)[*nabi + k].nname[0] = 'M';
257 else if (hbr->hack[j].tp == 11 && k >= 2)
259 /* This is a water lone pair, not a hydrogen */
260 /* Ugly hardcoded name hack */
261 srenew((*abi)[*nabi + k].nname, 4);
262 (*abi)[*nabi + k].nname[0] = 'L';
263 (*abi)[*nabi + k].nname[1] = 'P';
264 (*abi)[*nabi + k].nname[2] = '1' + k - 2;
265 (*abi)[*nabi + k].nname[3] = '\0';
267 else if (hbr->hack[j].nr > 1)
269 /* adding more than one atom, number them */
270 l = strlen((*abi)[*nabi + k].nname);
271 srenew((*abi)[*nabi + k].nname, l+2);
272 (*abi)[*nabi + k].nname[l] = '1' + k;
273 (*abi)[*nabi + k].nname[l+1] = '\0';
276 (*nabi) += hbr->hack[j].nr;
278 /* add hacks to atoms we've just added */
279 if (hbr->hack[j].tp > 0 || hbr->hack[j].oname == NULL)
281 for (k = 0; k < hbr->hack[j].nr; k++)
283 expand_hackblocks_one(hbr, (*abi)[*nabi-hbr->hack[j].nr+k].nname,
291 static void expand_hackblocks(t_atoms *pdba, t_hackblock hb[],
292 int nab[], t_hack *ab[],
293 int nterpairs, int *rN, int *rC)
298 for (i = 0; i < pdba->nr; i++)
301 for (j = 0; j < nterpairs && !bN; j++)
303 bN = pdba->atom[i].resind == rN[j];
306 for (j = 0; j < nterpairs && !bC; j++)
308 bC = pdba->atom[i].resind == rC[j];
311 /* add hacks to this atom */
312 expand_hackblocks_one(&hb[pdba->atom[i].resind], *pdba->atomname[i],
313 &nab[i], &ab[i], bN, bC);
317 fprintf(debug, "\n");
321 static int check_atoms_present(t_atoms *pdba, int nab[], t_hack *ab[])
323 int i, j, k, rnr, nadd;
326 for (i = 0; i < pdba->nr; i++)
328 rnr = pdba->atom[i].resind;
329 for (j = 0; j < nab[i]; j++)
331 if (ab[i][j].oname == NULL)
334 if (ab[i][j].nname == NULL)
336 gmx_incons("ab[i][j].nname not allocated");
338 /* check if the atom is already present */
339 k = pdbasearch_atom(ab[i][j].nname, rnr, pdba, "check", TRUE);
342 /* We found the added atom. */
343 ab[i][j].bAlreadyPresent = TRUE;
346 fprintf(debug, "Atom '%s' in residue '%s' %d is already present\n",
348 *pdba->resinfo[rnr].name, pdba->resinfo[rnr].nr);
353 ab[i][j].bAlreadyPresent = FALSE;
354 /* count how many atoms we'll add */
358 else if (ab[i][j].nname == NULL)
369 static void calc_all_pos(t_atoms *pdba, rvec x[], int nab[], t_hack *ab[],
370 gmx_bool bCheckMissing)
372 int i, j, ii, jj, m, ia, d, rnr, l = 0;
374 rvec xa[4]; /* control atoms for calc_h_pos */
375 rvec xh[MAXH]; /* hydrogen positions from calc_h_pos */
380 for (i = 0; i < pdba->nr; i++)
382 rnr = pdba->atom[i].resind;
383 for (j = 0; j < nab[i]; j += ab[i][j].nr)
385 /* check if we're adding: */
386 if (ab[i][j].oname == NULL && ab[i][j].tp > 0)
389 for (m = 0; (m < ab[i][j].nctl && bFoundAll); m++)
391 ia = pdbasearch_atom(ab[i][j].a[m], rnr, pdba,
392 bCheckMissing ? "atom" : "check",
396 /* not found in original atoms, might still be in t_hack (ab) */
397 hacksearch_atom(&ii, &jj, ab[i][j].a[m], nab, ab, rnr, pdba);
400 copy_rvec(ab[ii][jj].newx, xa[m]);
407 gmx_fatal(FARGS, "Atom %s not found in residue %s %d"
409 " while adding hydrogens",
411 *pdba->resinfo[rnr].name,
412 pdba->resinfo[rnr].nr,
413 *pdba->resinfo[rnr].rtp);
419 copy_rvec(x[ia], xa[m]);
424 for (m = 0; (m < MAXH); m++)
426 for (d = 0; d < DIM; d++)
438 calc_h_pos(ab[i][j].tp, xa, xh, &l);
439 for (m = 0; m < ab[i][j].nr; m++)
441 copy_rvec(xh[m], ab[i][j+m].newx);
442 ab[i][j+m].bXSet = TRUE;
450 static void free_ab(int natoms, int *nab, t_hack **ab)
454 for (i = 0; i < natoms; i++)
456 free_t_hack(nab[i], &ab[i]);
462 static int add_h_low(t_atoms **pdbaptr, rvec *xptr[],
463 int nah, t_hackblock ah[],
464 int nterpairs, t_hackblock **ntdb, t_hackblock **ctdb,
465 int *rN, int *rC, gmx_bool bCheckMissing,
466 int **nabptr, t_hack ***abptr,
467 gmx_bool bUpdate_pdba, gmx_bool bKeep_old_pdba)
469 t_atoms *newpdba = NULL, *pdba = NULL;
471 int i, newi, j, natoms, nalreadypresent;
478 /* set flags for adding hydrogens (according to hdb) */
484 /* the first time these will be pointers to NULL, but we will
485 return in them the completed arrays, which we will get back
492 fprintf(debug, "pointer to ab found\n");
502 /* WOW, everything was already figured out */
503 bUpdate_pdba = FALSE;
506 fprintf(debug, "pointer to non-null ab found\n");
511 /* We'll have to do all the hard work */
513 /* first get all the hackblocks for each residue: */
514 hb = get_hackblocks(pdba, nah, ah, nterpairs, ntdb, ctdb, rN, rC);
517 dump_hb(debug, pdba->nres, hb);
520 /* expand the hackblocks to atom level */
523 expand_hackblocks(pdba, hb, nab, ab, nterpairs, rN, rC);
524 free_t_hackblock(pdba->nres, &hb);
529 fprintf(debug, "before calc_all_pos\n");
530 dump_ab(debug, natoms, nab, ab, TRUE);
533 /* Now calc the positions */
534 calc_all_pos(pdba, *xptr, nab, ab, bCheckMissing);
538 fprintf(debug, "after calc_all_pos\n");
539 dump_ab(debug, natoms, nab, ab, TRUE);
544 /* we don't have to add atoms that are already present in pdba,
545 so we will remove them from the ab (t_hack) */
546 nadd = check_atoms_present(pdba, nab, ab);
549 fprintf(debug, "removed add hacks that were already in pdba:\n");
550 dump_ab(debug, natoms, nab, ab, TRUE);
551 fprintf(debug, "will be adding %d atoms\n", nadd);
554 /* Copy old atoms, making space for new ones */
556 init_t_atoms(newpdba, natoms+nadd, FALSE);
557 newpdba->nres = pdba->nres;
558 sfree(newpdba->resinfo);
559 newpdba->resinfo = pdba->resinfo;
567 fprintf(debug, "snew xn for %d old + %d new atoms %d total)\n",
568 natoms, nadd, natoms+nadd);
573 /* There is nothing to do: return now */
576 free_ab(natoms, nab, ab);
582 snew(xn, natoms+nadd);
584 for (i = 0; (i < natoms); i++)
586 /* check if this atom wasn't scheduled for deletion */
587 if (nab[i] == 0 || (ab[i][0].nname != NULL) )
589 if (newi >= natoms+nadd)
591 /*gmx_fatal(FARGS,"Not enough space for adding atoms");*/
593 srenew(xn, natoms+nadd);
596 srenew(newpdba->atom, natoms+nadd);
597 srenew(newpdba->atomname, natoms+nadd);
603 fprintf(debug, "(%3d) %3d %4s %4s%3d %3d",
604 i+1, newi+1, *pdba->atomname[i],
605 *pdba->resinfo[pdba->atom[i].resind].name,
606 pdba->resinfo[pdba->atom[i].resind].nr, nab[i]);
610 copy_atom(pdba, i, newpdba, newi);
612 copy_rvec((*xptr)[i], xn[newi]);
613 /* process the hacks for this atom */
615 for (j = 0; j < nab[i]; j++)
617 if (ab[i][j].oname == NULL) /* add */
620 if (newi >= natoms+nadd)
622 /* gmx_fatal(FARGS,"Not enough space for adding atoms");*/
624 srenew(xn, natoms+nadd);
627 srenew(newpdba->atom, natoms+nadd);
628 srenew(newpdba->atomname, natoms+nadd);
634 newpdba->atom[newi].resind = pdba->atom[i].resind;
638 fprintf(debug, " + %d", newi+1);
641 if (ab[i][j].nname != NULL &&
642 (ab[i][j].oname == NULL ||
643 strcmp(ab[i][j].oname, *newpdba->atomname[newi]) == 0))
646 if (ab[i][j].oname == NULL && ab[i][j].bAlreadyPresent)
648 /* This atom is already present, copy it from the input. */
652 copy_atom(pdba, i+nalreadypresent, newpdba, newi);
654 copy_rvec((*xptr)[i+nalreadypresent], xn[newi]);
662 fprintf(debug, "Replacing %d '%s' with (old name '%s') %s\n",
664 (newpdba->atomname[newi] && *newpdba->atomname[newi]) ? *newpdba->atomname[newi] : "",
665 ab[i][j].oname ? ab[i][j].oname : "",
668 snew(newpdba->atomname[newi], 1);
669 *newpdba->atomname[newi] = gmx_strdup(ab[i][j].nname);
670 if (ab[i][j].oname != NULL && ab[i][j].atom) /* replace */
671 { /* newpdba->atom[newi].m = ab[i][j].atom->m; */
672 /* newpdba->atom[newi].q = ab[i][j].atom->q; */
673 /* newpdba->atom[newi].type = ab[i][j].atom->type; */
678 copy_rvec(ab[i][j].newx, xn[newi]);
681 if (bUpdate_pdba && debug)
683 fprintf(debug, " %s %g %g", *newpdba->atomname[newi],
684 newpdba->atom[newi].m, newpdba->atom[newi].q);
689 i += nalreadypresent;
692 fprintf(debug, "\n");
709 free_ab(natoms, nab, ab);
716 for (i = 0; i < natoms; i++)
718 /* Do not free the atomname string itself, it might be in symtab */
719 /* sfree(*(pdba->atomname[i])); */
720 /* sfree(pdba->atomname[i]); */
722 sfree(pdba->atomname);
724 sfree(pdba->pdbinfo);
736 void deprotonate(t_atoms *atoms, rvec *x)
741 for (i = 0; i < atoms->nr; i++)
743 if ( (*atoms->atomname[i])[0] != 'H')
745 atoms->atomname[j] = atoms->atomname[i];
746 atoms->atom[j] = atoms->atom[i];
747 copy_rvec(x[i], x[j]);
754 int add_h(t_atoms **pdbaptr, rvec *xptr[],
755 int nah, t_hackblock ah[],
756 int nterpairs, t_hackblock **ntdb, t_hackblock **ctdb,
757 int *rN, int *rC, gmx_bool bAllowMissing,
758 int **nabptr, t_hack ***abptr,
759 gmx_bool bUpdate_pdba, gmx_bool bKeep_old_pdba)
761 int nold, nnew, niter;
763 /* Here we loop to be able to add atoms to added atoms.
764 * We should not check for missing atoms here.
771 nnew = add_h_low(pdbaptr, xptr, nah, ah, nterpairs, ntdb, ctdb, rN, rC, FALSE,
772 nabptr, abptr, bUpdate_pdba, bKeep_old_pdba);
776 gmx_fatal(FARGS, "More than 100 iterations of add_h. Maybe you are trying to replace an added atom (this is not supported)?");
783 /* Call add_h_low once more, now only for the missing atoms check */
784 add_h_low(pdbaptr, xptr, nah, ah, nterpairs, ntdb, ctdb, rN, rC, TRUE,
785 nabptr, abptr, bUpdate_pdba, bKeep_old_pdba);
791 int protonate(t_atoms **atomsptr, rvec **xptr, t_protonate *protdata)
795 gmx_bool bUpdate_pdba, bKeep_old_pdba;
796 int nntdb, nctdb, nt, ct;
800 if (!protdata->bInit)
804 fprintf(debug, "protonate: Initializing protdata\n");
807 /* set forcefield to use: */
808 strcpy(protdata->FF, "oplsaa.ff");
810 /* get the databases: */
811 protdata->nah = read_h_db(protdata->FF, &protdata->ah);
812 open_symtab(&protdata->tab);
813 protdata->atype = read_atype(protdata->FF, &protdata->tab);
814 nntdb = read_ter_db(protdata->FF, 'n', &protdata->ntdb, protdata->atype);
817 gmx_fatal(FARGS, "no N-terminus database");
819 nctdb = read_ter_db(protdata->FF, 'c', &protdata->ctdb, protdata->atype);
822 gmx_fatal(FARGS, "no C-terminus database");
825 /* set terminus types: -NH3+ (different for Proline) and -COO- */
827 snew(protdata->sel_ntdb, NTERPAIRS);
828 snew(protdata->sel_ctdb, NTERPAIRS);
830 if (nntdb >= 4 && nctdb >= 2)
832 /* Yuk, yuk, hardcoded default termini selections !!! */
833 if (strncmp(*atoms->resinfo[atoms->atom[atoms->nr-1].resind].name, "PRO", 3) == 0)
848 protdata->sel_ntdb[0] = &(protdata->ntdb[nt]);
849 protdata->sel_ctdb[0] = &(protdata->ctdb[ct]);
851 /* set terminal residue numbers: */
852 snew(protdata->rN, NTERPAIRS);
853 snew(protdata->rC, NTERPAIRS);
855 protdata->rC[0] = atoms->atom[atoms->nr-1].resind;
857 /* keep unprotonated topology: */
858 protdata->upatoms = atoms;
859 /* we don't have these yet: */
860 protdata->patoms = NULL;
862 bKeep_old_pdba = TRUE;
864 /* clear hackblocks: */
865 protdata->nab = NULL;
868 /* set flag to show we're initialized: */
869 protdata->bInit = TRUE;
875 fprintf(debug, "protonate: using available protdata\n");
877 /* add_h will need the unprotonated topology again: */
878 atoms = protdata->upatoms;
879 bUpdate_pdba = FALSE;
880 bKeep_old_pdba = FALSE;
884 nadd = add_h(&atoms, xptr, protdata->nah, protdata->ah,
885 NTERPAIRS, protdata->sel_ntdb, protdata->sel_ctdb,
886 protdata->rN, protdata->rC, TRUE,
887 &protdata->nab, &protdata->ab, bUpdate_pdba, bKeep_old_pdba);
888 if (!protdata->patoms)
890 /* store protonated topology */
891 protdata->patoms = atoms;
893 *atomsptr = protdata->patoms;
896 fprintf(debug, "natoms: %d -> %d (nadd=%d)\n",
897 protdata->upatoms->nr, protdata->patoms->nr, nadd);