change gmx_bool to bool in gmxpreprocess
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / genhydro.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "genhydro.h"
40
41 #include <string.h>
42 #include <time.h>
43
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"
58
59 static void copy_atom(t_atoms *atoms1, int a1, t_atoms *atoms2, int a2)
60 {
61     atoms2->atom[a2] = atoms1->atom[a1];
62     snew(atoms2->atomname[a2], 1);
63     *atoms2->atomname[a2] = gmx_strdup(*atoms1->atomname[a1]);
64 }
65
66 static int pdbasearch_atom(const char *name, int resind, t_atoms *pdba,
67                            const char *searchtype, bool bAllowMissing)
68 {
69     int  i;
70
71     for (i = 0; (i < pdba->nr) && (pdba->atom[i].resind != resind); i++)
72     {
73         ;
74     }
75
76     return search_atom(name, i, pdba,
77                        searchtype, bAllowMissing);
78 }
79
80 static void hacksearch_atom(int *ii, int *jj, char *name,
81                             const int nab[], t_hack *ab[],
82                             int resind, t_atoms *pdba)
83 {
84     int  i, j;
85
86     *ii = -1;
87     if (name[0] == '-')
88     {
89         name++;
90         resind--;
91     }
92     for (i = 0; (i < pdba->nr) && (pdba->atom[i].resind != resind); i++)
93     {
94         ;
95     }
96     for (; (i < pdba->nr) && (pdba->atom[i].resind == resind) && (*ii < 0); i++)
97     {
98         for (j = 0; (j < nab[i]) && (*ii < 0); j++)
99         {
100             if (ab[i][j].nname && strcmp(name, ab[i][j].nname) == 0)
101             {
102                 *ii = i;
103                 *jj = j;
104             }
105         }
106     }
107
108 }
109
110 static void dump_ab(FILE *out, int natom, const int nab[], t_hack *ab[], bool bHeader)
111 {
112     int i, j;
113
114 #define SS(s) (s) ? (s) : "-"
115     /* dump ab */
116     if (bHeader)
117     {
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");
121     }
122     for (i = 0; i < natom; i++)
123     {
124         for (j = 0; j < nab[i]; j++)
125         {
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),
128                     ab[i][j].tp,
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]);
133         }
134     }
135 #undef SS
136 }
137
138 static t_hackblock *get_hackblocks(t_atoms *pdba, int nah, t_hackblock ah[],
139                                    int nterpairs,
140                                    t_hackblock **ntdb, t_hackblock **ctdb,
141                                    const int *rN, const int *rC)
142 {
143     int          i, rnr;
144     t_hackblock *hb, *ahptr;
145
146     /* make space */
147     snew(hb, pdba->nres);
148     /* first the termini */
149     for (i = 0; i < nterpairs; i++)
150     {
151         if (ntdb[i] != nullptr)
152         {
153             copy_t_hackblock(ntdb[i], &hb[rN[i]]);
154         }
155         if (ctdb[i] != nullptr)
156         {
157             merge_t_hackblock(ctdb[i], &hb[rC[i]]);
158         }
159     }
160     /* then the whole hdb */
161     for (rnr = 0; rnr < pdba->nres; rnr++)
162     {
163         ahptr = search_h_db(nah, ah, *pdba->resinfo[rnr].rtp);
164         if (ahptr)
165         {
166             if (hb[rnr].name == nullptr)
167             {
168                 hb[rnr].name = gmx_strdup(ahptr->name);
169             }
170             merge_hacks(ahptr, &hb[rnr]);
171         }
172     }
173     return hb;
174 }
175
176 static void expand_hackblocks_one(t_hackblock *hbr, char *atomname,
177                                   int *nabi, t_hack **abi, bool bN, bool bC)
178 {
179     int      j, k, l;
180     bool     bIgnore;
181
182     /* we'll recursively add atoms to atoms */
183     for (j = 0; j < hbr->nhack; j++)
184     {
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) */
188         bIgnore = FALSE;
189         if (bN) /* N-terminus: ignore '-' */
190         {
191             for (k = 0; k < 4 && hbr->hack[j].a[k] && !bIgnore; k++)
192             {
193                 bIgnore = hbr->hack[j].a[k][0] == '-';
194             }
195         }
196         if (bC) /* C-terminus: ignore '+' */
197         {
198             for (k = 0; k < 4 && hbr->hack[j].a[k] && !bIgnore; k++)
199             {
200                 bIgnore = hbr->hack[j].a[k][0] == '+';
201             }
202         }
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 */
206         if (debug)
207         {
208             fprintf(debug, " %s", hbr->hack[j].oname ? hbr->hack[j].oname : hbr->hack[j].ai());
209         }
210
211         if (!bIgnore &&
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) ) )
216         {
217             /* now expand all hacks for this atom */
218             if (debug)
219             {
220                 fprintf(debug, " +%dh", hbr->hack[j].nr);
221             }
222             srenew(*abi, *nabi + hbr->hack[j].nr);
223             for (k = 0; k < hbr->hack[j].nr; k++)
224             {
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)
230                 {
231                     if ( (*abi)[*nabi + k].oname == nullptr)
232                     {
233                         (*abi)[*nabi + k].nname    = gmx_strdup(atomname);
234                         (*abi)[*nabi + k].nname[0] = 'H';
235                     }
236                 }
237                 else
238                 {
239                     if (gmx_debug_at)
240                     {
241                         fprintf(debug, "Hack '%s' %d, replacing nname '%s' with '%s' (old name '%s')\n",
242                                 atomname, j,
243                                 (*abi)[*nabi + k].nname, hbr->hack[j].nname,
244                                 (*abi)[*nabi + k].oname ? (*abi)[*nabi + k].oname : "");
245                     }
246                     sfree((*abi)[*nabi + k].nname);
247                     (*abi)[*nabi + k].nname = gmx_strdup(hbr->hack[j].nname);
248                 }
249
250                 if (hbr->hack[j].tp == 10 && k == 2)
251                 {
252                     /* This is a water virtual site, not a hydrogen */
253                     /* Ugly hardcoded name hack */
254                     (*abi)[*nabi + k].nname[0] = 'M';
255                 }
256                 else if (hbr->hack[j].tp == 11 && k >= 2)
257                 {
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';
265                 }
266                 else if (hbr->hack[j].nr > 1)
267                 {
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';
273                 }
274             }
275             (*nabi) += hbr->hack[j].nr;
276
277             /* add hacks to atoms we've just added */
278             if (hbr->hack[j].tp > 0 || hbr->hack[j].oname == nullptr)
279             {
280                 for (k = 0; k < hbr->hack[j].nr; k++)
281                 {
282                     expand_hackblocks_one(hbr, (*abi)[*nabi-hbr->hack[j].nr+k].nname,
283                                           nabi, abi, bN, bC);
284                 }
285             }
286         }
287     }
288 }
289
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)
293 {
294     int      i, j;
295     bool     bN, bC;
296
297     for (i = 0; i < pdba->nr; i++)
298     {
299         bN = FALSE;
300         for (j = 0; j < nterpairs && !bN; j++)
301         {
302             bN = pdba->atom[i].resind == rN[j];
303         }
304         bC = FALSE;
305         for (j = 0; j < nterpairs && !bC; j++)
306         {
307             bC = pdba->atom[i].resind == rC[j];
308         }
309
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);
313     }
314     if (debug)
315     {
316         fprintf(debug, "\n");
317     }
318 }
319
320 static int check_atoms_present(t_atoms *pdba, const int nab[], t_hack *ab[])
321 {
322     int i, j, k, rnr, nadd;
323
324     nadd = 0;
325     for (i = 0; i < pdba->nr; i++)
326     {
327         rnr = pdba->atom[i].resind;
328         for (j = 0; j < nab[i]; j++)
329         {
330             if (ab[i][j].oname == nullptr)
331             {
332                 /* we're adding */
333                 if (ab[i][j].nname == nullptr)
334                 {
335                     gmx_incons("ab[i][j].nname not allocated");
336                 }
337                 /* check if the atom is already present */
338                 k = pdbasearch_atom(ab[i][j].nname, rnr, pdba, "check", TRUE);
339                 if (k != -1)
340                 {
341                     /* We found the added atom. */
342                     ab[i][j].bAlreadyPresent = TRUE;
343                     if (debug)
344                     {
345                         fprintf(debug, "Atom '%s' in residue '%s' %d is already present\n",
346                                 ab[i][j].nname,
347                                 *pdba->resinfo[rnr].name, pdba->resinfo[rnr].nr);
348                     }
349                 }
350                 else
351                 {
352                     ab[i][j].bAlreadyPresent = FALSE;
353                     /* count how many atoms we'll add */
354                     nadd++;
355                 }
356             }
357             else if (ab[i][j].nname == nullptr)
358             {
359                 /* we're deleting */
360                 nadd--;
361             }
362         }
363     }
364
365     return nadd;
366 }
367
368 static void calc_all_pos(t_atoms *pdba, rvec x[], int nab[], t_hack *ab[],
369                          bool bCheckMissing)
370 {
371     int      i, j, ii, jj, m, ia, d, rnr, l = 0;
372 #define MAXH 4
373     rvec     xa[4];    /* control atoms for calc_h_pos */
374     rvec     xh[MAXH]; /* hydrogen positions from calc_h_pos */
375     bool     bFoundAll;
376
377     jj = 0;
378
379     for (i = 0; i < pdba->nr; i++)
380     {
381         rnr   = pdba->atom[i].resind;
382         for (j = 0; j < nab[i]; j += ab[i][j].nr)
383         {
384             /* check if we're adding: */
385             if (ab[i][j].oname == nullptr && ab[i][j].tp > 0)
386             {
387                 bFoundAll = TRUE;
388                 for (m = 0; (m < ab[i][j].nctl && bFoundAll); m++)
389                 {
390                     ia = pdbasearch_atom(ab[i][j].a[m], rnr, pdba,
391                                          bCheckMissing ? "atom" : "check",
392                                          !bCheckMissing);
393                     if (ia < 0)
394                     {
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);
397                         if (ii >= 0)
398                         {
399                             copy_rvec(ab[ii][jj].newx, xa[m]);
400                         }
401                         else
402                         {
403                             bFoundAll = FALSE;
404                             if (bCheckMissing)
405                             {
406                                 gmx_fatal(FARGS, "Atom %s not found in residue %s %d"
407                                           ", rtp entry %s"
408                                           " while adding hydrogens",
409                                           ab[i][j].a[m],
410                                           *pdba->resinfo[rnr].name,
411                                           pdba->resinfo[rnr].nr,
412                                           *pdba->resinfo[rnr].rtp);
413                             }
414                         }
415                     }
416                     else
417                     {
418                         copy_rvec(x[ia], xa[m]);
419                     }
420                 }
421                 if (bFoundAll)
422                 {
423                     for (m = 0; (m < MAXH); m++)
424                     {
425                         for (d = 0; d < DIM; d++)
426                         {
427                             if (m < ab[i][j].nr)
428                             {
429                                 xh[m][d] = 0;
430                             }
431                             else
432                             {
433                                 xh[m][d] = NOTSET;
434                             }
435                         }
436                     }
437                     calc_h_pos(ab[i][j].tp, xa, xh, &l);
438                     for (m = 0; m < ab[i][j].nr; m++)
439                     {
440                         copy_rvec(xh[m], ab[i][j+m].newx);
441                         ab[i][j+m].bXSet = TRUE;
442                     }
443                 }
444             }
445         }
446     }
447 }
448
449 static void free_ab(int natoms, int *nab, t_hack **ab)
450 {
451     int i;
452
453     for (i = 0; i < natoms; i++)
454     {
455         free_t_hack(nab[i], &ab[i]);
456     }
457     sfree(nab);
458     sfree(ab);
459 }
460
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)
467 {
468     t_atoms        *newpdba = nullptr, *pdba = nullptr;
469     int             nadd;
470     int             i, newi, j, natoms, nalreadypresent;
471     int            *nab = nullptr;
472     t_hack        **ab  = nullptr;
473     t_hackblock    *hb;
474     rvec           *xn;
475     bool            bKeep_ab;
476
477     /* set flags for adding hydrogens (according to hdb) */
478     pdba   = *pdbaptr;
479     natoms = pdba->nr;
480
481     if (nabptr && abptr)
482     {
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
485            the second time */
486         nab      = *nabptr;
487         ab       = *abptr;
488         bKeep_ab = TRUE;
489         if (debug)
490         {
491             fprintf(debug, "pointer to ab found\n");
492         }
493     }
494     else
495     {
496         bKeep_ab = FALSE;
497     }
498
499     if (nab && ab)
500     {
501         /* WOW, everything was already figured out */
502         bUpdate_pdba = FALSE;
503         if (debug)
504         {
505             fprintf(debug, "pointer to non-null ab found\n");
506         }
507     }
508     else
509     {
510         /* We'll have to do all the hard work */
511         bUpdate_pdba = TRUE;
512         /* first get all the hackblocks for each residue: */
513         hb = get_hackblocks(pdba, nah, ah, nterpairs, ntdb, ctdb, rN, rC);
514         if (debug)
515         {
516             dump_hb(debug, pdba->nres, hb);
517         }
518
519         /* expand the hackblocks to atom level */
520         snew(nab, natoms);
521         snew(ab, natoms);
522         expand_hackblocks(pdba, hb, nab, ab, nterpairs, rN, rC);
523         free_t_hackblock(pdba->nres, &hb);
524     }
525
526     if (debug)
527     {
528         fprintf(debug, "before calc_all_pos\n");
529         dump_ab(debug, natoms, nab, ab, TRUE);
530     }
531
532     /* Now calc the positions */
533     calc_all_pos(pdba, *xptr, nab, ab, bCheckMissing);
534
535     if (debug)
536     {
537         fprintf(debug, "after calc_all_pos\n");
538         dump_ab(debug, natoms, nab, ab, TRUE);
539     }
540
541     if (bUpdate_pdba)
542     {
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);
546         if (debug)
547         {
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);
551         }
552
553         /* Copy old atoms, making space for new ones */
554         snew(newpdba, 1);
555         init_t_atoms(newpdba, natoms+nadd, FALSE);
556         newpdba->nres    = pdba->nres;
557         sfree(newpdba->resinfo);
558         newpdba->resinfo = pdba->resinfo;
559     }
560     else
561     {
562         nadd = 0;
563     }
564     if (debug)
565     {
566         fprintf(debug, "snew xn for %d old + %d new atoms %d total)\n",
567                 natoms, nadd, natoms+nadd);
568     }
569
570     if (nadd == 0)
571     {
572         /* There is nothing to do: return now */
573         if (!bKeep_ab)
574         {
575             free_ab(natoms, nab, ab);
576         }
577
578         return natoms;
579     }
580
581     snew(xn, natoms+nadd);
582     newi = 0;
583     for (i = 0; (i < natoms); i++)
584     {
585         /* check if this atom wasn't scheduled for deletion */
586         if (nab[i] == 0 || (ab[i][0].nname != nullptr) )
587         {
588             if (newi >= natoms+nadd)
589             {
590                 /*gmx_fatal(FARGS,"Not enough space for adding atoms");*/
591                 nadd += 10;
592                 srenew(xn, natoms+nadd);
593                 if (bUpdate_pdba)
594                 {
595                     srenew(newpdba->atom, natoms+nadd);
596                     srenew(newpdba->atomname, natoms+nadd);
597                 }
598             }
599             if (debug)
600             {
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]);
605             }
606             if (bUpdate_pdba)
607             {
608                 copy_atom(pdba, i, newpdba, newi);
609             }
610             copy_rvec((*xptr)[i], xn[newi]);
611             /* process the hacks for this atom */
612             nalreadypresent = 0;
613             for (j = 0; j < nab[i]; j++)
614             {
615                 if (ab[i][j].oname == nullptr) /* add */
616                 {
617                     newi++;
618                     if (newi >= natoms+nadd)
619                     {
620                         /* gmx_fatal(FARGS,"Not enough space for adding atoms");*/
621                         nadd += 10;
622                         srenew(xn, natoms+nadd);
623                         if (bUpdate_pdba)
624                         {
625                             srenew(newpdba->atom, natoms+nadd);
626                             srenew(newpdba->atomname, natoms+nadd);
627                         }
628                     }
629                     if (bUpdate_pdba)
630                     {
631                         newpdba->atom[newi].resind = pdba->atom[i].resind;
632                     }
633                     if (debug)
634                     {
635                         fprintf(debug, " + %d", newi+1);
636                     }
637                 }
638                 if (ab[i][j].nname != nullptr &&
639                     (ab[i][j].oname == nullptr ||
640                      strcmp(ab[i][j].oname, *newpdba->atomname[newi]) == 0))
641                 {
642                     /* add or replace */
643                     if (ab[i][j].oname == nullptr && ab[i][j].bAlreadyPresent)
644                     {
645                         /* This atom is already present, copy it from the input. */
646                         nalreadypresent++;
647                         if (bUpdate_pdba)
648                         {
649                             copy_atom(pdba, i+nalreadypresent, newpdba, newi);
650                         }
651                         copy_rvec((*xptr)[i+nalreadypresent], xn[newi]);
652                     }
653                     else
654                     {
655                         if (bUpdate_pdba)
656                         {
657                             if (gmx_debug_at)
658                             {
659                                 fprintf(debug, "Replacing %d '%s' with (old name '%s') %s\n",
660                                         newi,
661                                         (newpdba->atomname[newi] && *newpdba->atomname[newi]) ? *newpdba->atomname[newi] : "",
662                                         ab[i][j].oname ? ab[i][j].oname : "",
663                                         ab[i][j].nname);
664                             }
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; */
671                             }
672                         }
673                         if (ab[i][j].bXSet)
674                         {
675                             copy_rvec(ab[i][j].newx, xn[newi]);
676                         }
677                     }
678                     if (bUpdate_pdba && debug)
679                     {
680                         fprintf(debug, " %s %g %g", *newpdba->atomname[newi],
681                                 newpdba->atom[newi].m, newpdba->atom[newi].q);
682                     }
683                 }
684             }
685             newi++;
686             i += nalreadypresent;
687             if (debug)
688             {
689                 fprintf(debug, "\n");
690             }
691         }
692     }
693     if (bUpdate_pdba)
694     {
695         newpdba->nr = newi;
696     }
697
698     if (bKeep_ab)
699     {
700         *nabptr = nab;
701         *abptr  = ab;
702     }
703     else
704     {
705         /* Clean up */
706         free_ab(natoms, nab, ab);
707     }
708
709     if (bUpdate_pdba)
710     {
711         if (!bKeep_old_pdba)
712         {
713             for (i = 0; i < natoms; i++)
714             {
715                 /* Do not free the atomname string itself, it might be in symtab */
716                 /* sfree(*(pdba->atomname[i])); */
717                 /* sfree(pdba->atomname[i]); */
718             }
719             sfree(pdba->atomname);
720             sfree(pdba->atom);
721             sfree(pdba->pdbinfo);
722             sfree(pdba);
723         }
724         *pdbaptr = newpdba;
725     }
726
727     sfree(*xptr);
728     *xptr = xn;
729
730     return newi;
731 }
732
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)
739 {
740     int nold, nnew, niter;
741
742     /* Here we loop to be able to add atoms to added atoms.
743      * We should not check for missing atoms here.
744      */
745     niter = 0;
746     nnew  = 0;
747     do
748     {
749         nold = nnew;
750         nnew = add_h_low(pdbaptr, xptr, nah, ah, nterpairs, ntdb, ctdb, rN, rC, FALSE,
751                          nabptr, abptr, bUpdate_pdba, bKeep_old_pdba);
752         niter++;
753         if (niter > 100)
754         {
755             gmx_fatal(FARGS, "More than 100 iterations of add_h. Maybe you are trying to replace an added atom (this is not supported)?");
756         }
757     }
758     while (nnew > nold);
759
760     if (!bAllowMissing)
761     {
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);
765     }
766
767     return nnew;
768 }