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,2017,2018, 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/gmxlib/network.h"
46 #include "gromacs/gmxpreprocess/calch.h"
47 #include "gromacs/gmxpreprocess/h_db.h"
48 #include "gromacs/gmxpreprocess/notset.h"
49 #include "gromacs/gmxpreprocess/pgutil.h"
50 #include "gromacs/gmxpreprocess/resall.h"
51 #include "gromacs/gmxpreprocess/ter_db.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 int pdbasearch_atom(const char *name, int resind, t_atoms *pdba,
67 const char *searchtype, 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 const 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)
110 static void dump_ab(FILE *out, int natom, const int nab[], t_hack *ab[], bool bHeader)
114 #define SS(s) (s) ? (s) : "-"
118 fprintf(out, "ADDBLOCK (t_hack) natom=%d\n"
119 "%4s %2s %-4s %-4s %2s %-4s %-4s %-4s %-4s %1s %s\n",
120 natom, "atom", "nr", "old", "new", "tp", "ai", "aj", "ak", "al", "a", "x");
122 for (i = 0; i < natom; i++)
124 for (j = 0; j < nab[i]; j++)
126 fprintf(out, "%4d %2d %-4s %-4s %2d %-4s %-4s %-4s %-4s %s %g %g %g\n",
127 i+1, ab[i][j].nr, SS(ab[i][j].oname), SS(ab[i][j].nname),
129 SS(ab[i][j].ai()), SS(ab[i][j].aj()),
130 SS(ab[i][j].ak()), SS(ab[i][j].al()),
131 ab[i][j].atom ? "+" : "",
132 ab[i][j].newx[XX], ab[i][j].newx[YY], ab[i][j].newx[ZZ]);
138 static t_hackblock *get_hackblocks(t_atoms *pdba, int nah, t_hackblock ah[],
140 t_hackblock **ntdb, t_hackblock **ctdb,
141 const int *rN, const int *rC)
144 t_hackblock *hb, *ahptr;
147 snew(hb, pdba->nres);
148 /* first the termini */
149 for (i = 0; i < nterpairs; i++)
151 if (ntdb[i] != nullptr)
153 copy_t_hackblock(ntdb[i], &hb[rN[i]]);
155 if (ctdb[i] != nullptr)
157 merge_t_hackblock(ctdb[i], &hb[rC[i]]);
160 /* then the whole hdb */
161 for (rnr = 0; rnr < pdba->nres; rnr++)
163 ahptr = search_h_db(nah, ah, *pdba->resinfo[rnr].rtp);
166 if (hb[rnr].name == nullptr)
168 hb[rnr].name = gmx_strdup(ahptr->name);
170 merge_hacks(ahptr, &hb[rnr]);
176 static void expand_hackblocks_one(t_hackblock *hbr, char *atomname,
177 int *nabi, t_hack **abi, bool bN, bool bC)
182 /* we'll recursively add atoms to atoms */
183 for (j = 0; j < hbr->nhack; j++)
185 /* first check if we're in the N- or C-terminus, then we should ignore
186 all hacks involving atoms from resp. previous or next residue
187 (i.e. which name begins with '-' (N) or '+' (C) */
189 if (bN) /* N-terminus: ignore '-' */
191 for (k = 0; k < 4 && hbr->hack[j].a[k] && !bIgnore; k++)
193 bIgnore = hbr->hack[j].a[k][0] == '-';
196 if (bC) /* C-terminus: ignore '+' */
198 for (k = 0; k < 4 && hbr->hack[j].a[k] && !bIgnore; k++)
200 bIgnore = hbr->hack[j].a[k][0] == '+';
203 /* must be either hdb entry (tp>0) or add from tdb (oname==NULL)
204 and first control aton (AI) matches this atom or
205 delete/replace from tdb (oname!=NULL) and oname matches this atom */
208 fprintf(debug, " %s", hbr->hack[j].oname ? hbr->hack[j].oname : hbr->hack[j].ai());
212 ( ( ( hbr->hack[j].tp > 0 || hbr->hack[j].oname == nullptr ) &&
213 strcmp(atomname, hbr->hack[j].ai()) == 0 ) ||
214 ( hbr->hack[j].oname != nullptr &&
215 strcmp(atomname, hbr->hack[j].oname) == 0) ) )
217 /* now expand all hacks for this atom */
220 fprintf(debug, " +%dh", hbr->hack[j].nr);
222 srenew(*abi, *nabi + hbr->hack[j].nr);
223 for (k = 0; k < hbr->hack[j].nr; k++)
225 copy_t_hack(&hbr->hack[j], &(*abi)[*nabi + k]);
226 (*abi)[*nabi + k].bXSet = FALSE;
227 /* if we're adding (oname==NULL) and don't have a new name (nname)
228 yet, build it from atomname */
229 if ( (*abi)[*nabi + k].nname == nullptr)
231 if ( (*abi)[*nabi + k].oname == nullptr)
233 (*abi)[*nabi + k].nname = gmx_strdup(atomname);
234 (*abi)[*nabi + k].nname[0] = 'H';
241 fprintf(debug, "Hack '%s' %d, replacing nname '%s' with '%s' (old name '%s')\n",
243 (*abi)[*nabi + k].nname, hbr->hack[j].nname,
244 (*abi)[*nabi + k].oname ? (*abi)[*nabi + k].oname : "");
246 sfree((*abi)[*nabi + k].nname);
247 (*abi)[*nabi + k].nname = gmx_strdup(hbr->hack[j].nname);
250 if (hbr->hack[j].tp == 10 && k == 2)
252 /* This is a water virtual site, not a hydrogen */
253 /* Ugly hardcoded name hack */
254 (*abi)[*nabi + k].nname[0] = 'M';
256 else if (hbr->hack[j].tp == 11 && k >= 2)
258 /* This is a water lone pair, not a hydrogen */
259 /* Ugly hardcoded name hack */
260 srenew((*abi)[*nabi + k].nname, 4);
261 (*abi)[*nabi + k].nname[0] = 'L';
262 (*abi)[*nabi + k].nname[1] = 'P';
263 (*abi)[*nabi + k].nname[2] = '1' + k - 2;
264 (*abi)[*nabi + k].nname[3] = '\0';
266 else if (hbr->hack[j].nr > 1)
268 /* adding more than one atom, number them */
269 l = strlen((*abi)[*nabi + k].nname);
270 srenew((*abi)[*nabi + k].nname, l+2);
271 (*abi)[*nabi + k].nname[l] = '1' + k;
272 (*abi)[*nabi + k].nname[l+1] = '\0';
275 (*nabi) += hbr->hack[j].nr;
277 /* add hacks to atoms we've just added */
278 if (hbr->hack[j].tp > 0 || hbr->hack[j].oname == nullptr)
280 for (k = 0; k < hbr->hack[j].nr; k++)
282 expand_hackblocks_one(hbr, (*abi)[*nabi-hbr->hack[j].nr+k].nname,
290 static void expand_hackblocks(t_atoms *pdba, t_hackblock hb[],
291 int nab[], t_hack *ab[],
292 int nterpairs, const int *rN, const int *rC)
297 for (i = 0; i < pdba->nr; i++)
300 for (j = 0; j < nterpairs && !bN; j++)
302 bN = pdba->atom[i].resind == rN[j];
305 for (j = 0; j < nterpairs && !bC; j++)
307 bC = pdba->atom[i].resind == rC[j];
310 /* add hacks to this atom */
311 expand_hackblocks_one(&hb[pdba->atom[i].resind], *pdba->atomname[i],
312 &nab[i], &ab[i], bN, bC);
316 fprintf(debug, "\n");
320 static int check_atoms_present(t_atoms *pdba, const int nab[], t_hack *ab[])
322 int i, j, k, rnr, nadd;
325 for (i = 0; i < pdba->nr; i++)
327 rnr = pdba->atom[i].resind;
328 for (j = 0; j < nab[i]; j++)
330 if (ab[i][j].oname == nullptr)
333 if (ab[i][j].nname == nullptr)
335 gmx_incons("ab[i][j].nname not allocated");
337 /* check if the atom is already present */
338 k = pdbasearch_atom(ab[i][j].nname, rnr, pdba, "check", TRUE);
341 /* We found the added atom. */
342 ab[i][j].bAlreadyPresent = TRUE;
345 fprintf(debug, "Atom '%s' in residue '%s' %d is already present\n",
347 *pdba->resinfo[rnr].name, pdba->resinfo[rnr].nr);
352 ab[i][j].bAlreadyPresent = FALSE;
353 /* count how many atoms we'll add */
357 else if (ab[i][j].nname == nullptr)
368 static void calc_all_pos(t_atoms *pdba, rvec x[], int nab[], t_hack *ab[],
371 int i, j, ii, jj, m, ia, d, rnr, l = 0;
373 rvec xa[4]; /* control atoms for calc_h_pos */
374 rvec xh[MAXH]; /* hydrogen positions from calc_h_pos */
379 for (i = 0; i < pdba->nr; i++)
381 rnr = pdba->atom[i].resind;
382 for (j = 0; j < nab[i]; j += ab[i][j].nr)
384 /* check if we're adding: */
385 if (ab[i][j].oname == nullptr && ab[i][j].tp > 0)
388 for (m = 0; (m < ab[i][j].nctl && bFoundAll); m++)
390 ia = pdbasearch_atom(ab[i][j].a[m], rnr, pdba,
391 bCheckMissing ? "atom" : "check",
395 /* not found in original atoms, might still be in t_hack (ab) */
396 hacksearch_atom(&ii, &jj, ab[i][j].a[m], nab, ab, rnr, pdba);
399 copy_rvec(ab[ii][jj].newx, xa[m]);
406 gmx_fatal(FARGS, "Atom %s not found in residue %s %d"
408 " while adding hydrogens",
410 *pdba->resinfo[rnr].name,
411 pdba->resinfo[rnr].nr,
412 *pdba->resinfo[rnr].rtp);
418 copy_rvec(x[ia], xa[m]);
423 for (m = 0; (m < MAXH); m++)
425 for (d = 0; d < DIM; d++)
437 calc_h_pos(ab[i][j].tp, xa, xh, &l);
438 for (m = 0; m < ab[i][j].nr; m++)
440 copy_rvec(xh[m], ab[i][j+m].newx);
441 ab[i][j+m].bXSet = TRUE;
449 static void free_ab(int natoms, int *nab, t_hack **ab)
453 for (i = 0; i < natoms; i++)
455 free_t_hack(nab[i], &ab[i]);
461 static int add_h_low(t_atoms **pdbaptr, rvec *xptr[],
462 int nah, t_hackblock ah[],
463 int nterpairs, t_hackblock **ntdb, t_hackblock **ctdb,
464 int *rN, int *rC, bool bCheckMissing,
465 int **nabptr, t_hack ***abptr,
466 bool bUpdate_pdba, bool bKeep_old_pdba)
468 t_atoms *newpdba = nullptr, *pdba = nullptr;
470 int i, newi, j, natoms, nalreadypresent;
472 t_hack **ab = nullptr;
477 /* set flags for adding hydrogens (according to hdb) */
483 /* the first time these will be pointers to NULL, but we will
484 return in them the completed arrays, which we will get back
491 fprintf(debug, "pointer to ab found\n");
501 /* WOW, everything was already figured out */
502 bUpdate_pdba = FALSE;
505 fprintf(debug, "pointer to non-null ab found\n");
510 /* We'll have to do all the hard work */
512 /* first get all the hackblocks for each residue: */
513 hb = get_hackblocks(pdba, nah, ah, nterpairs, ntdb, ctdb, rN, rC);
516 dump_hb(debug, pdba->nres, hb);
519 /* expand the hackblocks to atom level */
522 expand_hackblocks(pdba, hb, nab, ab, nterpairs, rN, rC);
523 free_t_hackblock(pdba->nres, &hb);
528 fprintf(debug, "before calc_all_pos\n");
529 dump_ab(debug, natoms, nab, ab, TRUE);
532 /* Now calc the positions */
533 calc_all_pos(pdba, *xptr, nab, ab, bCheckMissing);
537 fprintf(debug, "after calc_all_pos\n");
538 dump_ab(debug, natoms, nab, ab, TRUE);
543 /* we don't have to add atoms that are already present in pdba,
544 so we will remove them from the ab (t_hack) */
545 nadd = check_atoms_present(pdba, nab, ab);
548 fprintf(debug, "removed add hacks that were already in pdba:\n");
549 dump_ab(debug, natoms, nab, ab, TRUE);
550 fprintf(debug, "will be adding %d atoms\n", nadd);
553 /* Copy old atoms, making space for new ones */
555 init_t_atoms(newpdba, natoms+nadd, FALSE);
556 newpdba->nres = pdba->nres;
557 sfree(newpdba->resinfo);
558 newpdba->resinfo = pdba->resinfo;
566 fprintf(debug, "snew xn for %d old + %d new atoms %d total)\n",
567 natoms, nadd, natoms+nadd);
572 /* There is nothing to do: return now */
575 free_ab(natoms, nab, ab);
581 snew(xn, natoms+nadd);
583 for (i = 0; (i < natoms); i++)
585 /* check if this atom wasn't scheduled for deletion */
586 if (nab[i] == 0 || (ab[i][0].nname != nullptr) )
588 if (newi >= natoms+nadd)
590 /*gmx_fatal(FARGS,"Not enough space for adding atoms");*/
592 srenew(xn, natoms+nadd);
595 srenew(newpdba->atom, natoms+nadd);
596 srenew(newpdba->atomname, natoms+nadd);
601 fprintf(debug, "(%3d) %3d %4s %4s%3d %3d",
602 i+1, newi+1, *pdba->atomname[i],
603 *pdba->resinfo[pdba->atom[i].resind].name,
604 pdba->resinfo[pdba->atom[i].resind].nr, nab[i]);
608 copy_atom(pdba, i, newpdba, newi);
610 copy_rvec((*xptr)[i], xn[newi]);
611 /* process the hacks for this atom */
613 for (j = 0; j < nab[i]; j++)
615 if (ab[i][j].oname == nullptr) /* add */
618 if (newi >= natoms+nadd)
620 /* gmx_fatal(FARGS,"Not enough space for adding atoms");*/
622 srenew(xn, natoms+nadd);
625 srenew(newpdba->atom, natoms+nadd);
626 srenew(newpdba->atomname, natoms+nadd);
631 newpdba->atom[newi].resind = pdba->atom[i].resind;
635 fprintf(debug, " + %d", newi+1);
638 if (ab[i][j].nname != nullptr &&
639 (ab[i][j].oname == nullptr ||
640 strcmp(ab[i][j].oname, *newpdba->atomname[newi]) == 0))
643 if (ab[i][j].oname == nullptr && ab[i][j].bAlreadyPresent)
645 /* This atom is already present, copy it from the input. */
649 copy_atom(pdba, i+nalreadypresent, newpdba, newi);
651 copy_rvec((*xptr)[i+nalreadypresent], xn[newi]);
659 fprintf(debug, "Replacing %d '%s' with (old name '%s') %s\n",
661 (newpdba->atomname[newi] && *newpdba->atomname[newi]) ? *newpdba->atomname[newi] : "",
662 ab[i][j].oname ? ab[i][j].oname : "",
665 snew(newpdba->atomname[newi], 1);
666 *newpdba->atomname[newi] = gmx_strdup(ab[i][j].nname);
667 if (ab[i][j].oname != nullptr && ab[i][j].atom) /* replace */
668 { /* newpdba->atom[newi].m = ab[i][j].atom->m; */
669 /* newpdba->atom[newi].q = ab[i][j].atom->q; */
670 /* newpdba->atom[newi].type = ab[i][j].atom->type; */
675 copy_rvec(ab[i][j].newx, xn[newi]);
678 if (bUpdate_pdba && debug)
680 fprintf(debug, " %s %g %g", *newpdba->atomname[newi],
681 newpdba->atom[newi].m, newpdba->atom[newi].q);
686 i += nalreadypresent;
689 fprintf(debug, "\n");
706 free_ab(natoms, nab, ab);
713 for (i = 0; i < natoms; i++)
715 /* Do not free the atomname string itself, it might be in symtab */
716 /* sfree(*(pdba->atomname[i])); */
717 /* sfree(pdba->atomname[i]); */
719 sfree(pdba->atomname);
721 sfree(pdba->pdbinfo);
733 int add_h(t_atoms **pdbaptr, rvec *xptr[],
734 int nah, t_hackblock ah[],
735 int nterpairs, t_hackblock **ntdb, t_hackblock **ctdb,
736 int *rN, int *rC, bool bAllowMissing,
737 int **nabptr, t_hack ***abptr,
738 bool bUpdate_pdba, bool bKeep_old_pdba)
740 int nold, nnew, niter;
742 /* Here we loop to be able to add atoms to added atoms.
743 * We should not check for missing atoms here.
750 nnew = add_h_low(pdbaptr, xptr, nah, ah, nterpairs, ntdb, ctdb, rN, rC, FALSE,
751 nabptr, abptr, bUpdate_pdba, bKeep_old_pdba);
755 gmx_fatal(FARGS, "More than 100 iterations of add_h. Maybe you are trying to replace an added atom (this is not supported)?");
762 /* Call add_h_low once more, now only for the missing atoms check */
763 add_h_low(pdbaptr, xptr, nah, ah, nterpairs, ntdb, ctdb, rN, rC, TRUE,
764 nabptr, abptr, bUpdate_pdba, bKeep_old_pdba);