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/notset.h"
48 #include "gromacs/gmxpreprocess/pgutil.h"
49 #include "gromacs/gmxpreprocess/resall.h"
50 #include "gromacs/gmxpreprocess/ter_db.h"
51 #include "gromacs/legacyheaders/network.h"
52 #include "gromacs/legacyheaders/typedefs.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/topology/symtab.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
60 static void copy_atom(t_atoms
*atoms1
, int a1
, t_atoms
*atoms2
, int a2
)
62 atoms2
->atom
[a2
] = atoms1
->atom
[a1
];
63 snew(atoms2
->atomname
[a2
], 1);
64 *atoms2
->atomname
[a2
] = gmx_strdup(*atoms1
->atomname
[a1
]);
67 static atom_id
pdbasearch_atom(const char *name
, int resind
, t_atoms
*pdba
,
68 const char *searchtype
, gmx_bool bAllowMissing
)
72 for (i
= 0; (i
< pdba
->nr
) && (pdba
->atom
[i
].resind
!= resind
); i
++)
77 return search_atom(name
, i
, pdba
,
78 searchtype
, bAllowMissing
);
81 static void hacksearch_atom(int *ii
, int *jj
, char *name
,
82 int nab
[], t_hack
*ab
[],
83 int resind
, t_atoms
*pdba
)
93 for (i
= 0; (i
< pdba
->nr
) && (pdba
->atom
[i
].resind
!= resind
); i
++)
97 for (; (i
< pdba
->nr
) && (pdba
->atom
[i
].resind
== resind
) && (*ii
< 0); i
++)
99 for (j
= 0; (j
< nab
[i
]) && (*ii
< 0); j
++)
101 if (ab
[i
][j
].nname
&& strcmp(name
, ab
[i
][j
].nname
) == 0)
112 void dump_ab(FILE *out
, int natom
, int nab
[], t_hack
*ab
[], gmx_bool bHeader
)
116 #define SS(s) (s) ? (s) : "-"
120 fprintf(out
, "ADDBLOCK (t_hack) natom=%d\n"
121 "%4s %2s %-4s %-4s %2s %-4s %-4s %-4s %-4s %1s %s\n",
122 natom
, "atom", "nr", "old", "new", "tp", "ai", "aj", "ak", "al", "a", "x");
124 for (i
= 0; i
< natom
; i
++)
126 for (j
= 0; j
< nab
[i
]; j
++)
128 fprintf(out
, "%4d %2d %-4s %-4s %2d %-4s %-4s %-4s %-4s %s %g %g %g\n",
129 i
+1, ab
[i
][j
].nr
, SS(ab
[i
][j
].oname
), SS(ab
[i
][j
].nname
),
131 SS(ab
[i
][j
].ai()), SS(ab
[i
][j
].aj()),
132 SS(ab
[i
][j
].ak()), SS(ab
[i
][j
].al()),
133 ab
[i
][j
].atom
? "+" : "",
134 ab
[i
][j
].newx
[XX
], ab
[i
][j
].newx
[YY
], ab
[i
][j
].newx
[ZZ
]);
140 static t_hackblock
*get_hackblocks(t_atoms
*pdba
, int nah
, t_hackblock ah
[],
142 t_hackblock
**ntdb
, t_hackblock
**ctdb
,
146 t_hackblock
*hb
, *ahptr
;
149 snew(hb
, pdba
->nres
);
150 /* first the termini */
151 for (i
= 0; i
< nterpairs
; i
++)
155 copy_t_hackblock(ntdb
[i
], &hb
[rN
[i
]]);
159 merge_t_hackblock(ctdb
[i
], &hb
[rC
[i
]]);
162 /* then the whole hdb */
163 for (rnr
= 0; rnr
< pdba
->nres
; rnr
++)
165 ahptr
= search_h_db(nah
, ah
, *pdba
->resinfo
[rnr
].rtp
);
168 if (hb
[rnr
].name
== NULL
)
170 hb
[rnr
].name
= gmx_strdup(ahptr
->name
);
172 merge_hacks(ahptr
, &hb
[rnr
]);
178 static void expand_hackblocks_one(t_hackblock
*hbr
, char *atomname
,
179 int *nabi
, t_hack
**abi
, gmx_bool bN
, gmx_bool bC
)
184 /* we'll recursively add atoms to atoms */
185 for (j
= 0; j
< hbr
->nhack
; j
++)
187 /* first check if we're in the N- or C-terminus, then we should ignore
188 all hacks involving atoms from resp. previous or next residue
189 (i.e. which name begins with '-' (N) or '+' (C) */
191 if (bN
) /* N-terminus: ignore '-' */
193 for (k
= 0; k
< 4 && hbr
->hack
[j
].a
[k
] && !bIgnore
; k
++)
195 bIgnore
= hbr
->hack
[j
].a
[k
][0] == '-';
198 if (bC
) /* C-terminus: ignore '+' */
200 for (k
= 0; k
< 4 && hbr
->hack
[j
].a
[k
] && !bIgnore
; k
++)
202 bIgnore
= hbr
->hack
[j
].a
[k
][0] == '+';
205 /* must be either hdb entry (tp>0) or add from tdb (oname==NULL)
206 and first control aton (AI) matches this atom or
207 delete/replace from tdb (oname!=NULL) and oname matches this atom */
210 fprintf(debug
, " %s", hbr
->hack
[j
].oname
? hbr
->hack
[j
].oname
: hbr
->hack
[j
].ai());
214 ( ( ( hbr
->hack
[j
].tp
> 0 || hbr
->hack
[j
].oname
== NULL
) &&
215 strcmp(atomname
, hbr
->hack
[j
].ai()) == 0 ) ||
216 ( hbr
->hack
[j
].oname
!= NULL
&&
217 strcmp(atomname
, hbr
->hack
[j
].oname
) == 0) ) )
219 /* now expand all hacks for this atom */
222 fprintf(debug
, " +%dh", hbr
->hack
[j
].nr
);
224 srenew(*abi
, *nabi
+ hbr
->hack
[j
].nr
);
225 for (k
= 0; k
< hbr
->hack
[j
].nr
; k
++)
227 copy_t_hack(&hbr
->hack
[j
], &(*abi
)[*nabi
+ k
]);
228 (*abi
)[*nabi
+ k
].bXSet
= FALSE
;
229 /* if we're adding (oname==NULL) and don't have a new name (nname)
230 yet, build it from atomname */
231 if ( (*abi
)[*nabi
+ k
].nname
== NULL
)
233 if ( (*abi
)[*nabi
+ k
].oname
== NULL
)
235 (*abi
)[*nabi
+ k
].nname
= gmx_strdup(atomname
);
236 (*abi
)[*nabi
+ k
].nname
[0] = 'H';
243 fprintf(debug
, "Hack '%s' %d, replacing nname '%s' with '%s' (old name '%s')\n",
245 (*abi
)[*nabi
+ k
].nname
, hbr
->hack
[j
].nname
,
246 (*abi
)[*nabi
+ k
].oname
? (*abi
)[*nabi
+ k
].oname
: "");
248 sfree((*abi
)[*nabi
+ k
].nname
);
249 (*abi
)[*nabi
+ k
].nname
= gmx_strdup(hbr
->hack
[j
].nname
);
252 if (hbr
->hack
[j
].tp
== 10 && k
== 2)
254 /* This is a water virtual site, not a hydrogen */
255 /* Ugly hardcoded name hack */
256 (*abi
)[*nabi
+ k
].nname
[0] = 'M';
258 else if (hbr
->hack
[j
].tp
== 11 && k
>= 2)
260 /* This is a water lone pair, not a hydrogen */
261 /* Ugly hardcoded name hack */
262 srenew((*abi
)[*nabi
+ k
].nname
, 4);
263 (*abi
)[*nabi
+ k
].nname
[0] = 'L';
264 (*abi
)[*nabi
+ k
].nname
[1] = 'P';
265 (*abi
)[*nabi
+ k
].nname
[2] = '1' + k
- 2;
266 (*abi
)[*nabi
+ k
].nname
[3] = '\0';
268 else if (hbr
->hack
[j
].nr
> 1)
270 /* adding more than one atom, number them */
271 l
= strlen((*abi
)[*nabi
+ k
].nname
);
272 srenew((*abi
)[*nabi
+ k
].nname
, l
+2);
273 (*abi
)[*nabi
+ k
].nname
[l
] = '1' + k
;
274 (*abi
)[*nabi
+ k
].nname
[l
+1] = '\0';
277 (*nabi
) += hbr
->hack
[j
].nr
;
279 /* add hacks to atoms we've just added */
280 if (hbr
->hack
[j
].tp
> 0 || hbr
->hack
[j
].oname
== NULL
)
282 for (k
= 0; k
< hbr
->hack
[j
].nr
; k
++)
284 expand_hackblocks_one(hbr
, (*abi
)[*nabi
-hbr
->hack
[j
].nr
+k
].nname
,
292 static void expand_hackblocks(t_atoms
*pdba
, t_hackblock hb
[],
293 int nab
[], t_hack
*ab
[],
294 int nterpairs
, int *rN
, int *rC
)
299 for (i
= 0; i
< pdba
->nr
; i
++)
302 for (j
= 0; j
< nterpairs
&& !bN
; j
++)
304 bN
= pdba
->atom
[i
].resind
== rN
[j
];
307 for (j
= 0; j
< nterpairs
&& !bC
; j
++)
309 bC
= pdba
->atom
[i
].resind
== rC
[j
];
312 /* add hacks to this atom */
313 expand_hackblocks_one(&hb
[pdba
->atom
[i
].resind
], *pdba
->atomname
[i
],
314 &nab
[i
], &ab
[i
], bN
, bC
);
318 fprintf(debug
, "\n");
322 static int check_atoms_present(t_atoms
*pdba
, int nab
[], t_hack
*ab
[])
324 int i
, j
, k
, rnr
, nadd
;
327 for (i
= 0; i
< pdba
->nr
; i
++)
329 rnr
= pdba
->atom
[i
].resind
;
330 for (j
= 0; j
< nab
[i
]; j
++)
332 if (ab
[i
][j
].oname
== NULL
)
335 if (ab
[i
][j
].nname
== NULL
)
337 gmx_incons("ab[i][j].nname not allocated");
339 /* check if the atom is already present */
340 k
= pdbasearch_atom(ab
[i
][j
].nname
, rnr
, pdba
, "check", TRUE
);
343 /* We found the added atom. */
344 ab
[i
][j
].bAlreadyPresent
= TRUE
;
347 fprintf(debug
, "Atom '%s' in residue '%s' %d is already present\n",
349 *pdba
->resinfo
[rnr
].name
, pdba
->resinfo
[rnr
].nr
);
354 ab
[i
][j
].bAlreadyPresent
= FALSE
;
355 /* count how many atoms we'll add */
359 else if (ab
[i
][j
].nname
== NULL
)
370 static void calc_all_pos(t_atoms
*pdba
, rvec x
[], int nab
[], t_hack
*ab
[],
371 gmx_bool bCheckMissing
)
373 int i
, j
, ii
, jj
, m
, ia
, d
, rnr
, l
= 0;
375 rvec xa
[4]; /* control atoms for calc_h_pos */
376 rvec xh
[MAXH
]; /* hydrogen positions from calc_h_pos */
381 for (i
= 0; i
< pdba
->nr
; i
++)
383 rnr
= pdba
->atom
[i
].resind
;
384 for (j
= 0; j
< nab
[i
]; j
+= ab
[i
][j
].nr
)
386 /* check if we're adding: */
387 if (ab
[i
][j
].oname
== NULL
&& ab
[i
][j
].tp
> 0)
390 for (m
= 0; (m
< ab
[i
][j
].nctl
&& bFoundAll
); m
++)
392 ia
= pdbasearch_atom(ab
[i
][j
].a
[m
], rnr
, pdba
,
393 bCheckMissing
? "atom" : "check",
397 /* not found in original atoms, might still be in t_hack (ab) */
398 hacksearch_atom(&ii
, &jj
, ab
[i
][j
].a
[m
], nab
, ab
, rnr
, pdba
);
401 copy_rvec(ab
[ii
][jj
].newx
, xa
[m
]);
408 gmx_fatal(FARGS
, "Atom %s not found in residue %s %d"
410 " while adding hydrogens",
412 *pdba
->resinfo
[rnr
].name
,
413 pdba
->resinfo
[rnr
].nr
,
414 *pdba
->resinfo
[rnr
].rtp
);
420 copy_rvec(x
[ia
], xa
[m
]);
425 for (m
= 0; (m
< MAXH
); m
++)
427 for (d
= 0; d
< DIM
; d
++)
439 calc_h_pos(ab
[i
][j
].tp
, xa
, xh
, &l
);
440 for (m
= 0; m
< ab
[i
][j
].nr
; m
++)
442 copy_rvec(xh
[m
], ab
[i
][j
+m
].newx
);
443 ab
[i
][j
+m
].bXSet
= TRUE
;
451 static void free_ab(int natoms
, int *nab
, t_hack
**ab
)
455 for (i
= 0; i
< natoms
; i
++)
457 free_t_hack(nab
[i
], &ab
[i
]);
463 static int add_h_low(t_atoms
**pdbaptr
, rvec
*xptr
[],
464 int nah
, t_hackblock ah
[],
465 int nterpairs
, t_hackblock
**ntdb
, t_hackblock
**ctdb
,
466 int *rN
, int *rC
, gmx_bool bCheckMissing
,
467 int **nabptr
, t_hack
***abptr
,
468 gmx_bool bUpdate_pdba
, gmx_bool bKeep_old_pdba
)
470 t_atoms
*newpdba
= NULL
, *pdba
= NULL
;
472 int i
, newi
, j
, natoms
, nalreadypresent
;
479 /* set flags for adding hydrogens (according to hdb) */
485 /* the first time these will be pointers to NULL, but we will
486 return in them the completed arrays, which we will get back
493 fprintf(debug
, "pointer to ab found\n");
503 /* WOW, everything was already figured out */
504 bUpdate_pdba
= FALSE
;
507 fprintf(debug
, "pointer to non-null ab found\n");
512 /* We'll have to do all the hard work */
514 /* first get all the hackblocks for each residue: */
515 hb
= get_hackblocks(pdba
, nah
, ah
, nterpairs
, ntdb
, ctdb
, rN
, rC
);
518 dump_hb(debug
, pdba
->nres
, hb
);
521 /* expand the hackblocks to atom level */
524 expand_hackblocks(pdba
, hb
, nab
, ab
, nterpairs
, rN
, rC
);
525 free_t_hackblock(pdba
->nres
, &hb
);
530 fprintf(debug
, "before calc_all_pos\n");
531 dump_ab(debug
, natoms
, nab
, ab
, TRUE
);
534 /* Now calc the positions */
535 calc_all_pos(pdba
, *xptr
, nab
, ab
, bCheckMissing
);
539 fprintf(debug
, "after calc_all_pos\n");
540 dump_ab(debug
, natoms
, nab
, ab
, TRUE
);
545 /* we don't have to add atoms that are already present in pdba,
546 so we will remove them from the ab (t_hack) */
547 nadd
= check_atoms_present(pdba
, nab
, ab
);
550 fprintf(debug
, "removed add hacks that were already in pdba:\n");
551 dump_ab(debug
, natoms
, nab
, ab
, TRUE
);
552 fprintf(debug
, "will be adding %d atoms\n", nadd
);
555 /* Copy old atoms, making space for new ones */
557 init_t_atoms(newpdba
, natoms
+nadd
, FALSE
);
558 newpdba
->nres
= pdba
->nres
;
559 sfree(newpdba
->resinfo
);
560 newpdba
->resinfo
= pdba
->resinfo
;
568 fprintf(debug
, "snew xn for %d old + %d new atoms %d total)\n",
569 natoms
, nadd
, natoms
+nadd
);
574 /* There is nothing to do: return now */
577 free_ab(natoms
, nab
, ab
);
583 snew(xn
, natoms
+nadd
);
585 for (i
= 0; (i
< natoms
); i
++)
587 /* check if this atom wasn't scheduled for deletion */
588 if (nab
[i
] == 0 || (ab
[i
][0].nname
!= NULL
) )
590 if (newi
>= natoms
+nadd
)
592 /*gmx_fatal(FARGS,"Not enough space for adding atoms");*/
594 srenew(xn
, natoms
+nadd
);
597 srenew(newpdba
->atom
, natoms
+nadd
);
598 srenew(newpdba
->atomname
, natoms
+nadd
);
604 fprintf(debug
, "(%3d) %3d %4s %4s%3d %3d",
605 i
+1, newi
+1, *pdba
->atomname
[i
],
606 *pdba
->resinfo
[pdba
->atom
[i
].resind
].name
,
607 pdba
->resinfo
[pdba
->atom
[i
].resind
].nr
, nab
[i
]);
611 copy_atom(pdba
, i
, newpdba
, newi
);
613 copy_rvec((*xptr
)[i
], xn
[newi
]);
614 /* process the hacks for this atom */
616 for (j
= 0; j
< nab
[i
]; j
++)
618 if (ab
[i
][j
].oname
== NULL
) /* add */
621 if (newi
>= natoms
+nadd
)
623 /* gmx_fatal(FARGS,"Not enough space for adding atoms");*/
625 srenew(xn
, natoms
+nadd
);
628 srenew(newpdba
->atom
, natoms
+nadd
);
629 srenew(newpdba
->atomname
, natoms
+nadd
);
635 newpdba
->atom
[newi
].resind
= pdba
->atom
[i
].resind
;
639 fprintf(debug
, " + %d", newi
+1);
642 if (ab
[i
][j
].nname
!= NULL
&&
643 (ab
[i
][j
].oname
== NULL
||
644 strcmp(ab
[i
][j
].oname
, *newpdba
->atomname
[newi
]) == 0))
647 if (ab
[i
][j
].oname
== NULL
&& ab
[i
][j
].bAlreadyPresent
)
649 /* This atom is already present, copy it from the input. */
653 copy_atom(pdba
, i
+nalreadypresent
, newpdba
, newi
);
655 copy_rvec((*xptr
)[i
+nalreadypresent
], xn
[newi
]);
663 fprintf(debug
, "Replacing %d '%s' with (old name '%s') %s\n",
665 (newpdba
->atomname
[newi
] && *newpdba
->atomname
[newi
]) ? *newpdba
->atomname
[newi
] : "",
666 ab
[i
][j
].oname
? ab
[i
][j
].oname
: "",
669 snew(newpdba
->atomname
[newi
], 1);
670 *newpdba
->atomname
[newi
] = gmx_strdup(ab
[i
][j
].nname
);
671 if (ab
[i
][j
].oname
!= NULL
&& ab
[i
][j
].atom
) /* replace */
672 { /* newpdba->atom[newi].m = ab[i][j].atom->m; */
673 /* newpdba->atom[newi].q = ab[i][j].atom->q; */
674 /* newpdba->atom[newi].type = ab[i][j].atom->type; */
679 copy_rvec(ab
[i
][j
].newx
, xn
[newi
]);
682 if (bUpdate_pdba
&& debug
)
684 fprintf(debug
, " %s %g %g", *newpdba
->atomname
[newi
],
685 newpdba
->atom
[newi
].m
, newpdba
->atom
[newi
].q
);
690 i
+= nalreadypresent
;
693 fprintf(debug
, "\n");
710 free_ab(natoms
, nab
, ab
);
717 for (i
= 0; i
< natoms
; i
++)
719 /* Do not free the atomname string itself, it might be in symtab */
720 /* sfree(*(pdba->atomname[i])); */
721 /* sfree(pdba->atomname[i]); */
723 sfree(pdba
->atomname
);
725 sfree(pdba
->pdbinfo
);
737 int add_h(t_atoms
**pdbaptr
, rvec
*xptr
[],
738 int nah
, t_hackblock ah
[],
739 int nterpairs
, t_hackblock
**ntdb
, t_hackblock
**ctdb
,
740 int *rN
, int *rC
, gmx_bool bAllowMissing
,
741 int **nabptr
, t_hack
***abptr
,
742 gmx_bool bUpdate_pdba
, gmx_bool bKeep_old_pdba
)
744 int nold
, nnew
, niter
;
746 /* Here we loop to be able to add atoms to added atoms.
747 * We should not check for missing atoms here.
754 nnew
= add_h_low(pdbaptr
, xptr
, nah
, ah
, nterpairs
, ntdb
, ctdb
, rN
, rC
, FALSE
,
755 nabptr
, abptr
, bUpdate_pdba
, bKeep_old_pdba
);
759 gmx_fatal(FARGS
, "More than 100 iterations of add_h. Maybe you are trying to replace an added atom (this is not supported)?");
766 /* Call add_h_low once more, now only for the missing atoms check */
767 add_h_low(pdbaptr
, xptr
, nah
, ah
, nterpairs
, ntdb
, ctdb
, rN
, rC
, TRUE
,
768 nabptr
, abptr
, bUpdate_pdba
, bKeep_old_pdba
);