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