OpenMM: added combination rule check, disabled restraint check for now as it's buggy
[gromacs/qmmm-gamess-us.git] / src / kernel / genhydro.c
blob5250d20aad620a0f282700e0387170e23f630f95
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
40 #include <time.h>
41 #include <ctype.h>
42 #include "sysstuff.h"
43 #include "typedefs.h"
44 #include "smalloc.h"
45 #include "copyrite.h"
46 #include "string2.h"
47 #include "confio.h"
48 #include "symtab.h"
49 #include "vec.h"
50 #include "statutil.h"
51 #include "futil.h"
52 #include "gmx_fatal.h"
53 #include "physics.h"
54 #include "calch.h"
55 #include "genhydro.h"
56 #include "h_db.h"
57 #include "ter_db.h"
58 #include "resall.h"
59 #include "pgutil.h"
60 #include "network.h"
62 static void copy_atom(t_atoms *atoms1,int a1,t_atoms *atoms2,int a2)
64 atoms2->atom[a2] = atoms1->atom[a1];
65 snew(atoms2->atomname[a2],1);
66 *atoms2->atomname[a2]=strdup(*atoms1->atomname[a1]);
69 static atom_id pdbasearch_atom(const char *name,int resind,t_atoms *pdba,
70 const char *searchtype,bool bAllowMissing)
72 int i;
74 for(i=0; (i<pdba->nr) && (pdba->atom[i].resind != resind); i++)
77 return search_atom(name,i,pdba->nr,pdba->atom,pdba->atomname,
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)
85 int i,j;
87 *ii=-1;
88 if (name[0] == '-') {
89 name++;
90 resind--;
92 for(i=0; (i<pdba->nr) && (pdba->atom[i].resind != resind); i++)
94 for( ; (i<pdba->nr) && (pdba->atom[i].resind == resind) && (*ii<0); i++)
95 for(j=0; (j < nab[i]) && (*ii<0); j++)
96 if (ab[i][j].nname && strcmp(name,ab[i][j].nname) == 0) {
97 *ii=i;
98 *jj=j;
101 return;
104 void dump_ab(FILE *out,int natom,int nab[], t_hack *ab[], bool bHeader)
106 int i,j;
108 #define SS(s) (s)?(s):"-"
109 /* dump ab */
110 if (bHeader)
111 fprintf(out,"ADDBLOCK (t_hack) natom=%d\n"
112 "%4s %2s %-4s %-4s %2s %-4s %-4s %-4s %-4s %1s %s\n",
113 natom,"atom","nr","old","new","tp","ai","aj","ak","al","a","x");
114 for(i=0; i<natom; i++)
115 for(j=0; j<nab[i]; j++)
116 fprintf(out,"%4d %2d %-4s %-4s %2d %-4s %-4s %-4s %-4s %s %g %g %g\n",
117 i+1, ab[i][j].nr, SS(ab[i][j].oname), SS(ab[i][j].nname),
118 ab[i][j].tp,
119 SS(ab[i][j].AI), SS(ab[i][j].AJ),
120 SS(ab[i][j].AK), SS(ab[i][j].AL),
121 ab[i][j].atom?"+":"",
122 ab[i][j].newx[XX], ab[i][j].newx[YY], ab[i][j].newx[ZZ]);
123 #undef SS
126 static t_hackblock *get_hackblocks(t_atoms *pdba, int nah, t_hackblock ah[],
127 int nterpairs,
128 t_hackblock **ntdb, t_hackblock **ctdb,
129 int *rN, int *rC)
131 int i,rnr;
132 t_hackblock *hb,*ahptr;
134 /* make space */
135 snew(hb,pdba->nres);
136 /* first the termini */
137 for(i=0; i<nterpairs; i++) {
138 if (ntdb[i] != NULL) {
139 copy_t_hackblock(ntdb[i], &hb[rN[i]]);
141 if (ctdb[i] != NULL) {
142 merge_t_hackblock(ctdb[i], &hb[rC[i]]);
145 /* then the whole hdb */
146 for(rnr=0; rnr < pdba->nres; rnr++) {
147 ahptr=search_h_db(nah,ah,*pdba->resinfo[rnr].rtp);
148 if ( ahptr ) {
149 if (hb[rnr].name==NULL)
150 hb[rnr].name=strdup(ahptr->name);
151 merge_hacks(ahptr, &hb[rnr]);
154 return hb;
157 static const char Hnum[] = "123456";
159 static void expand_hackblocks_one(t_hackblock *hbr, char *atomname,
160 int *nabi, t_hack **abi, bool bN, bool bC)
162 int j, k, l, d;
163 bool bIgnore;
165 /* we'll recursively add atoms to atoms */
166 for(j=0; j < hbr->nhack; j++) {
167 /* first check if we're in the N- or C-terminus, then we should ignore
168 all hacks involving atoms from resp. previous or next residue
169 (i.e. which name begins with '-' (N) or '+' (C) */
170 bIgnore = FALSE;
171 if ( bN ) /* N-terminus: ignore '-' */
172 for(k=0; k<4 && hbr->hack[j].a[k] && !bIgnore; k++)
173 bIgnore = hbr->hack[j].a[k][0]=='-';
174 if ( bC ) /* C-terminus: ignore '+' */
175 for(k=0; k<4 && hbr->hack[j].a[k] && !bIgnore; k++)
176 bIgnore = hbr->hack[j].a[k][0]=='+';
177 /* must be either hdb entry (tp>0) or add from tdb (oname==NULL)
178 and first control aton (AI) matches this atom or
179 delete/replace from tdb (oname!=NULL) and oname matches this atom */
180 if (debug) fprintf(debug," %s",
181 hbr->hack[j].oname?hbr->hack[j].oname:hbr->hack[j].AI);
183 if ( !bIgnore &&
184 ( ( ( hbr->hack[j].tp > 0 || hbr->hack[j].oname==NULL ) &&
185 strcmp(atomname, hbr->hack[j].AI) == 0 ) ||
186 ( hbr->hack[j].oname!=NULL &&
187 strcmp(atomname, hbr->hack[j].oname) == 0) ) ) {
188 /* now expand all hacks for this atom */
189 if (debug) fprintf(debug," +%dh",hbr->hack[j].nr);
190 srenew(*abi,*nabi + hbr->hack[j].nr);
191 for(k=0; k < hbr->hack[j].nr; k++) {
192 copy_t_hack(&hbr->hack[j], &(*abi)[*nabi + k]);
193 (*abi)[*nabi + k].bXSet = FALSE;
194 /* if we're adding (oname==NULL) and don't have a new name (nname)
195 yet, build it from atomname */
196 if ( (*abi)[*nabi + k].nname==NULL ) {
197 if ( (*abi)[*nabi + k].oname==NULL ) {
198 (*abi)[*nabi + k].nname=strdup(atomname);
199 (*abi)[*nabi + k].nname[0]='H';
201 } else {
202 if (gmx_debug_at) {
203 fprintf(debug,"Hack '%s' %d, replacing nname '%s' with '%s' (old name '%s')\n",
204 atomname,j,
205 (*abi)[*nabi + k].nname,hbr->hack[j].nname,
206 (*abi)[*nabi + k].oname ? (*abi)[*nabi + k].oname : "");
208 sfree((*abi)[*nabi + k].nname);
209 (*abi)[*nabi + k].nname=strdup(hbr->hack[j].nname);
211 /* if adding more than one atom, number them */
212 if ( hbr->hack[j].nr > 1 ) {
213 l = strlen((*abi)[*nabi + k].nname);
214 srenew((*abi)[*nabi + k].nname, l+2);
215 (*abi)[*nabi + k].nname[l] = Hnum[k]; /* 1, 2, 3 .... */
216 (*abi)[*nabi + k].nname[l+1] = '\0';
219 (*nabi) += hbr->hack[j].nr;
221 /* add hacks to atoms we've just added */
222 if ( hbr->hack[j].tp > 0 || hbr->hack[j].oname==NULL )
223 for(k=0; k < hbr->hack[j].nr; k++)
224 expand_hackblocks_one(hbr, (*abi)[*nabi-hbr->hack[j].nr+k].nname,
225 nabi, abi, bN, bC);
230 static void expand_hackblocks(t_atoms *pdba, t_hackblock hb[],
231 int nab[], t_hack *ab[],
232 int nterpairs, int *rN, int *rC)
234 int i,j;
235 bool bN,bC;
237 for(i=0; i < pdba->nr; i++) {
238 bN = FALSE;
239 for(j=0; j<nterpairs && !bN; j++)
240 bN = pdba->atom[i].resind==rN[j];
241 bC = FALSE;
242 for(j=0; j<nterpairs && !bC; j++)
243 bC = pdba->atom[i].resind==rC[j];
245 /* add hacks to this atom */
246 expand_hackblocks_one(&hb[pdba->atom[i].resind], *pdba->atomname[i],
247 &nab[i], &ab[i], bN, bC);
249 if (debug) fprintf(debug,"\n");
252 static int check_atoms_present(t_atoms *pdba, int nab[], t_hack *ab[])
254 int i, j, k, d, rnr, nadd;
256 nadd=0;
257 for(i=0; i < pdba->nr; i++) {
258 rnr = pdba->atom[i].resind;
259 for(j=0; j<nab[i]; j++) {
260 if ( ab[i][j].oname==NULL ) {
261 /* we're adding */
262 if (ab[i][j].nname == NULL)
263 gmx_incons("ab[i][j].nname not allocated");
264 /* check if the atom is already present */
265 k = pdbasearch_atom(ab[i][j].nname, rnr, pdba, "check", TRUE);
266 if ( k != -1 ) {
267 /* We found the added atom. */
268 ab[i][j].bAlreadyPresent = TRUE;
269 if (debug) {
270 fprintf(debug,"Atom '%s' in residue '%s' %d is already present\n",
271 ab[i][j].nname,
272 *pdba->resinfo[rnr].name,pdba->resinfo[rnr].nr);
274 } else {
275 ab[i][j].bAlreadyPresent = FALSE;
276 /* count how many atoms we'll add */
277 nadd++;
279 } else if ( ab[i][j].nname==NULL ) {
280 /* we're deleting */
281 nadd--;
286 return nadd;
289 static void calc_all_pos(t_atoms *pdba, rvec x[], int nab[], t_hack *ab[],
290 bool bCheckMissing)
292 int i, j, ii, jj, m, ia, d, rnr,l=0;
293 #define MAXH 4
294 rvec xa[4]; /* control atoms for calc_h_pos */
295 rvec xh[MAXH]; /* hydrogen positions from calc_h_pos */
296 bool bFoundAll;
298 jj = 0;
300 for(i=0; i < pdba->nr; i++) {
301 rnr = pdba->atom[i].resind;
302 for(j=0; j < nab[i]; j+=ab[i][j].nr) {
303 /* check if we're adding: */
304 if (ab[i][j].oname==NULL && ab[i][j].tp > 0) {
305 bFoundAll = TRUE;
306 for(m=0; (m<ab[i][j].nctl && bFoundAll); m++) {
307 ia = pdbasearch_atom(ab[i][j].a[m], rnr, pdba,
308 bCheckMissing ? "atom" : "check",
309 !bCheckMissing);
310 if (ia < 0) {
311 /* not found in original atoms, might still be in t_hack (ab) */
312 hacksearch_atom(&ii, &jj, ab[i][j].a[m], nab, ab, rnr, pdba);
313 if (ii >= 0) {
314 copy_rvec(ab[ii][jj].newx, xa[m]);
315 } else {
316 bFoundAll = FALSE;
317 if (bCheckMissing) {
318 gmx_fatal(FARGS,"Atom %s not found in residue %s %d"
319 ", rtp entry %s"
320 " while adding hydrogens",
321 ab[i][j].a[m],
322 *pdba->resinfo[rnr].name,
323 pdba->resinfo[rnr].nr,
324 *pdba->resinfo[rnr].rtp);
327 } else
328 copy_rvec(x[ia], xa[m]);
330 if (bFoundAll) {
331 for(m=0; (m<MAXH); m++)
332 for(d=0; d<DIM; d++)
333 if (m<ab[i][j].nr)
334 xh[m][d] = 0;
335 else
336 xh[m][d] = NOTSET;
337 calc_h_pos(ab[i][j].tp, xa, xh, &l);
338 for(m=0; m<ab[i][j].nr; m++) {
339 copy_rvec(xh[m],ab[i][j+m].newx);
340 ab[i][j+m].bXSet = TRUE;
348 static void free_ab(int natoms,int *nab,t_hack **ab)
350 int i;
352 for(i=0; i<natoms; i++) {
353 free_t_hack(nab[i], &ab[i]);
355 sfree(nab);
356 sfree(ab);
359 static int add_h_low(t_atoms **pdbaptr, rvec *xptr[],
360 int nah, t_hackblock ah[],
361 int nterpairs, t_hackblock **ntdb, t_hackblock **ctdb,
362 int *rN, int *rC, bool bCheckMissing,
363 int **nabptr, t_hack ***abptr,
364 bool bUpdate_pdba, bool bKeep_old_pdba)
366 t_atoms *newpdba=NULL,*pdba=NULL;
367 int nadd;
368 int i,newi,j,d,natoms,nalreadypresent;
369 int *nab=NULL;
370 t_hack **ab=NULL;
371 t_hackblock *hb;
372 rvec *xn;
373 bool bKeep_ab;
375 /* set flags for adding hydrogens (according to hdb) */
376 pdba=*pdbaptr;
377 natoms=pdba->nr;
379 if (nabptr && abptr) {
380 /* the first time these will be pointers to NULL, but we will
381 return in them the completed arrays, which we will get back
382 the second time */
383 nab = *nabptr;
384 ab = *abptr;
385 bKeep_ab=TRUE;
386 if (debug) fprintf(debug,"pointer to ab found\n");
387 } else {
388 bKeep_ab=FALSE;
391 if (nab && ab) {
392 /* WOW, everything was already figured out */
393 bUpdate_pdba = FALSE;
394 if (debug) fprintf(debug,"pointer to non-null ab found\n");
395 } else {
396 /* We'll have to do all the hard work */
397 bUpdate_pdba = TRUE;
398 /* first get all the hackblocks for each residue: */
399 hb = get_hackblocks(pdba, nah, ah, nterpairs, ntdb, ctdb, rN, rC);
400 if (debug) dump_hb(debug, pdba->nres, hb);
402 /* expand the hackblocks to atom level */
403 snew(nab,natoms);
404 snew(ab,natoms);
405 expand_hackblocks(pdba, hb, nab, ab, nterpairs, rN, rC);
406 free_t_hackblock(pdba->nres, &hb);
409 if (debug) {
410 fprintf(debug,"before calc_all_pos\n");
411 dump_ab(debug, natoms, nab, ab, TRUE);
414 /* Now calc the positions */
415 calc_all_pos(pdba, *xptr, nab, ab, bCheckMissing);
417 if (debug) {
418 fprintf(debug,"after calc_all_pos\n");
419 dump_ab(debug, natoms, nab, ab, TRUE);
422 if (bUpdate_pdba) {
423 /* we don't have to add atoms that are already present in pdba,
424 so we will remove them from the ab (t_hack) */
425 nadd = check_atoms_present(pdba, nab, ab);
426 if (debug) {
427 fprintf(debug, "removed add hacks that were already in pdba:\n");
428 dump_ab(debug, natoms, nab, ab, TRUE);
429 fprintf(debug, "will be adding %d atoms\n",nadd);
432 /* Copy old atoms, making space for new ones */
433 snew(newpdba,1);
434 init_t_atoms(newpdba,natoms+nadd,FALSE);
435 newpdba->nres = pdba->nres;
436 sfree(newpdba->resinfo);
437 newpdba->resinfo = pdba->resinfo;
438 } else {
439 nadd = 0;
441 if (debug) fprintf(debug,"snew xn for %d old + %d new atoms %d total)\n",
442 natoms, nadd, natoms+nadd);
444 if (nadd == 0) {
445 /* There is nothing to do: return now */
446 if (!bKeep_ab) {
447 free_ab(natoms,nab,ab);
450 return natoms;
453 snew(xn,natoms+nadd);
454 newi=0;
455 for(i=0; (i<natoms); i++) {
456 /* check if this atom wasn't scheduled for deletion */
457 if ( nab[i]==0 || (ab[i][0].nname != NULL) ) {
458 if (newi >= natoms+nadd) {
459 /*gmx_fatal(FARGS,"Not enough space for adding atoms");*/
460 nadd+=10;
461 srenew(xn,natoms+nadd);
462 if (bUpdate_pdba) {
463 srenew(newpdba->atom,natoms+nadd);
464 srenew(newpdba->atomname,natoms+nadd);
466 debug_gmx();
468 if (debug) fprintf(debug,"(%3d) %3d %4s %4s%3d %3d",
469 i+1, newi+1, *pdba->atomname[i],
470 *pdba->resinfo[pdba->atom[i].resind].name,
471 pdba->resinfo[pdba->atom[i].resind].nr, nab[i]);
472 if (bUpdate_pdba)
473 copy_atom(pdba,i, newpdba,newi);
474 copy_rvec((*xptr)[i],xn[newi]);
475 /* process the hacks for this atom */
476 nalreadypresent = 0;
477 for(j=0; j<nab[i]; j++) {
478 if ( ab[i][j].oname==NULL ) { /* add */
479 newi++;
480 if (newi >= natoms+nadd) {
481 /* gmx_fatal(FARGS,"Not enough space for adding atoms");*/
482 nadd+=10;
483 srenew(xn,natoms+nadd);
484 if (bUpdate_pdba) {
485 srenew(newpdba->atom,natoms+nadd);
486 srenew(newpdba->atomname,natoms+nadd);
488 debug_gmx();
490 if (bUpdate_pdba) {
491 newpdba->atom[newi].resind=pdba->atom[i].resind;
493 if (debug) fprintf(debug," + %d",newi+1);
495 if (ab[i][j].nname != NULL &&
496 (ab[i][j].oname == NULL ||
497 strcmp(ab[i][j].oname,*newpdba->atomname[newi]) == 0)) {
498 /* add or replace */
499 if (ab[i][j].oname == NULL && ab[i][j].bAlreadyPresent) {
500 /* This atom is already present, copy it from the input. */
501 nalreadypresent++;
502 if (bUpdate_pdba) {
503 copy_atom(pdba,i+nalreadypresent, newpdba,newi);
505 copy_rvec((*xptr)[i+nalreadypresent],xn[newi]);
506 } else {
507 if (bUpdate_pdba) {
508 if (gmx_debug_at) {
509 fprintf(debug,"Replacing %d '%s' with (old name '%s') %s\n",
510 newi,
511 (newpdba->atomname[newi] && *newpdba->atomname[newi]) ? *newpdba->atomname[newi] : "",
512 ab[i][j].oname ? ab[i][j].oname : "",
513 ab[i][j].nname);
515 snew(newpdba->atomname[newi],1);
516 *newpdba->atomname[newi]=strdup(ab[i][j].nname);
517 if ( ab[i][j].oname!=NULL && ab[i][j].atom ) { /* replace */
518 /* newpdba->atom[newi].m = ab[i][j].atom->m; */
519 /* newpdba->atom[newi].q = ab[i][j].atom->q; */
520 /* newpdba->atom[newi].type = ab[i][j].atom->type; */
523 if (ab[i][j].bXSet)
524 copy_rvec(ab[i][j].newx, xn[newi]);
526 if (bUpdate_pdba && debug)
527 fprintf(debug," %s %g %g",*newpdba->atomname[newi],
528 newpdba->atom[newi].m,newpdba->atom[newi].q);
531 newi++;
532 i += nalreadypresent;
533 if (debug) fprintf(debug,"\n");
536 if (bUpdate_pdba)
537 newpdba->nr = newi;
539 if ( bKeep_ab ) {
540 *nabptr=nab;
541 *abptr=ab;
542 } else {
543 /* Clean up */
544 free_ab(natoms,nab,ab);
547 if ( bUpdate_pdba ) {
548 if ( !bKeep_old_pdba ) {
549 for(i=0; i < natoms; i++) {
550 /* Do not free the atomname string itself, it might be in symtab */
551 /* sfree(*(pdba->atomname[i])); */
552 /* sfree(pdba->atomname[i]); */
554 sfree(pdba->atomname);
555 sfree(pdba->atom);
556 sfree(pdba->pdbinfo);
557 sfree(pdba);
559 *pdbaptr=newpdba;
560 } else
561 nadd = newi-natoms;
563 sfree(*xptr);
564 *xptr=xn;
566 return newi;
569 void deprotonate(t_atoms *atoms,rvec *x)
571 int i,j;
573 j=0;
574 for(i=0; i<atoms->nr; i++) {
575 if ( (*atoms->atomname[i])[0] != 'H' ) {
576 atoms->atomname[j]=atoms->atomname[i];
577 atoms->atom[j]=atoms->atom[i];
578 copy_rvec(x[i],x[j]);
579 j++;
582 atoms->nr=j;
585 int add_h(t_atoms **pdbaptr, rvec *xptr[],
586 int nah, t_hackblock ah[],
587 int nterpairs, t_hackblock **ntdb, t_hackblock **ctdb,
588 int *rN, int *rC, bool bAllowMissing,
589 int **nabptr, t_hack ***abptr,
590 bool bUpdate_pdba, bool bKeep_old_pdba)
592 int nold,nnew,niter;
594 /* Here we loop to be able to add atoms to added atoms.
595 * We should not check for missing atoms here.
597 niter = 0;
598 nnew = 0;
599 do {
600 nold = nnew;
601 nnew = add_h_low(pdbaptr,xptr,nah,ah,nterpairs,ntdb,ctdb,rN,rC,FALSE,
602 nabptr,abptr,bUpdate_pdba,bKeep_old_pdba);
603 niter++;
604 if (niter > 100) {
605 gmx_fatal(FARGS,"More than 100 iterations of add_h. Maybe you are trying to replace an added atom (this is not supported)?");
607 } while (nnew > nold);
609 if (!bAllowMissing) {
610 /* Call add_h_low once more, now only for the missing atoms check */
611 add_h_low(pdbaptr,xptr,nah,ah,nterpairs,ntdb,ctdb,rN,rC,TRUE,
612 nabptr,abptr,bUpdate_pdba,bKeep_old_pdba);
615 return nnew;
618 int protonate(t_atoms **atomsptr,rvec **xptr,t_protonate *protdata)
620 #define NTERPAIRS 1
621 t_atoms *atoms;
622 bool bUpdate_pdba,bKeep_old_pdba;
623 int nntdb,nctdb,nt,ct;
624 int nadd;
626 atoms=NULL;
627 if ( !protdata->bInit ) {
628 if (debug) fprintf(debug,"protonate: Initializing protdata\n");
629 /* set forcefield to use: */
630 strcpy(protdata->FF,"ffgmx2");
631 /* get the databases: */
632 protdata->nah=read_h_db(protdata->FF,&protdata->ah);
633 open_symtab(&protdata->tab);
634 protdata->atype=read_atype(protdata->FF,&protdata->tab);
635 nntdb = read_ter_db(protdata->FF,'n',&protdata->ntdb,protdata->atype);
636 if (nntdb < 1)
637 gmx_fatal(FARGS,"no n-terminus db");
638 nctdb = read_ter_db(protdata->FF,'c',&protdata->ctdb,protdata->atype);
639 if (nctdb < 1)
640 gmx_fatal(FARGS,"no c-terminus db");
642 /* set terminus types: -NH3+ (different for Proline) and -COO- */
643 atoms=*atomsptr;
644 snew(protdata->sel_ntdb, NTERPAIRS);
645 snew(protdata->sel_ctdb, NTERPAIRS);
647 if (nntdb>=4 && nctdb>=2) {
648 /* Yuk, yuk, hardcoded default termini selections !!! */
649 if (strncmp(*atoms->resinfo[atoms->atom[atoms->nr-1].resind].name,"PRO",3)==0)
650 nt = 3;
651 else
652 nt = 1;
653 ct = 1;
654 } else {
655 nt = 0;
656 ct = 0;
658 protdata->sel_ntdb[0]=&(protdata->ntdb[nt]);
659 protdata->sel_ctdb[0]=&(protdata->ctdb[ct]);
661 /* set terminal residue numbers: */
662 snew(protdata->rN, NTERPAIRS);
663 snew(protdata->rC, NTERPAIRS);
664 protdata->rN[0]=0;
665 protdata->rC[0]=atoms->atom[atoms->nr-1].resind;
667 /* keep unprotonated topology: */
668 protdata->upatoms = atoms;
669 /* we don't have these yet: */
670 protdata->patoms = NULL;
671 bUpdate_pdba=TRUE;
672 bKeep_old_pdba=TRUE;
674 /* clear hackblocks: */
675 protdata->nab=NULL;
676 protdata->ab=NULL;
678 /* set flag to show we're initialized: */
679 protdata->bInit=TRUE;
680 } else {
681 if (debug) fprintf(debug,"protonate: using available protdata\n");
682 /* add_h will need the unprotonated topoloy again: */
683 atoms = protdata->upatoms;
684 bUpdate_pdba=FALSE;
685 bKeep_old_pdba=FALSE;
688 /* now protonate */
689 nadd = add_h(&atoms, xptr, protdata->nah, protdata->ah,
690 NTERPAIRS, protdata->sel_ntdb, protdata->sel_ctdb,
691 protdata->rN, protdata->rC, TRUE,
692 &protdata->nab, &protdata->ab, bUpdate_pdba, bKeep_old_pdba);
693 if ( ! protdata->patoms )
694 /* store protonated topology */
695 protdata->patoms = atoms;
696 *atomsptr = protdata->patoms;
697 if (debug) fprintf(debug,"natoms: %d -> %d (nadd=%d)\n",
698 protdata->upatoms->nr, protdata->patoms->nr, nadd);
699 return nadd;
700 #undef NTERPAIRS