Replaced NOTSET from typedefs.h by local solutions.
[gromacs.git] / src / gromacs / gmxpreprocess / gpp_atomtype.cpp
blob67bf8da643570549f1047a4bb100474e8363db7a
1 /*
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) 2011,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.
37 #include "gmxpre.h"
39 #include "gpp_atomtype.h"
41 #include <string.h>
43 #include <climits>
44 #include <cmath>
46 #include "gromacs/gmxpreprocess/notset.h"
47 #include "gromacs/gmxpreprocess/topdirs.h"
48 #include "gromacs/gmxpreprocess/toputil.h"
49 #include "gromacs/legacyheaders/txtdump.h"
50 #include "gromacs/legacyheaders/types/ifunc.h"
51 #include "gromacs/topology/symtab.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/smalloc.h"
55 typedef struct gpp_atomtype {
56 int nr; /* The number of atomtypes */
57 t_atom *atom; /* Array of atoms */
58 char ***atomname; /* Names of the atomtypes */
59 t_param *nb; /* Nonbonded force default params */
60 int *bondatomtype; /* The bond_atomtype for each atomtype */
61 real *radius; /* Radius for GBSA stuff */
62 real *vol; /* Effective volume for GBSA */
63 real *surftens; /* Surface tension with water, for GBSA */
64 real *gb_radius; /* Radius for Still model */
65 real *S_hct; /* Overlap factor for HCT model */
66 int *atomnumber; /* Atomic number, used for QM/MM */
67 } t_gpp_atomtype;
69 int get_atomtype_type(const char *str, gpp_atomtype_t ga)
71 int i;
73 /* Atom types are always case sensitive */
74 for (i = 0; (i < ga->nr); i++)
76 if (strcmp(str, *(ga->atomname[i])) == 0)
78 return i;
82 return NOTSET;
85 int get_atomtype_ntypes(gpp_atomtype_t ga)
87 return ga->nr;
90 char *get_atomtype_name(int nt, gpp_atomtype_t ga)
92 if ((nt < 0) || (nt >= ga->nr))
94 return NULL;
97 return *(ga->atomname[nt]);
100 real get_atomtype_massA(int nt, gpp_atomtype_t ga)
102 if ((nt < 0) || (nt >= ga->nr))
104 return NOTSET;
107 return ga->atom[nt].m;
110 real get_atomtype_massB(int nt, gpp_atomtype_t ga)
112 if ((nt < 0) || (nt >= ga->nr))
114 return NOTSET;
117 return ga->atom[nt].mB;
120 real get_atomtype_qA(int nt, gpp_atomtype_t ga)
122 if ((nt < 0) || (nt >= ga->nr))
124 return NOTSET;
127 return ga->atom[nt].q;
130 real get_atomtype_qB(int nt, gpp_atomtype_t ga)
132 if ((nt < 0) || (nt >= ga->nr))
134 return NOTSET;
137 return ga->atom[nt].qB;
140 int get_atomtype_ptype(int nt, gpp_atomtype_t ga)
142 if ((nt < 0) || (nt >= ga->nr))
144 return NOTSET;
147 return ga->atom[nt].ptype;
150 int get_atomtype_batype(int nt, gpp_atomtype_t ga)
152 if ((nt < 0) || (nt >= ga->nr))
154 return NOTSET;
157 return ga->bondatomtype[nt];
160 int get_atomtype_atomnumber(int nt, gpp_atomtype_t ga)
162 if ((nt < 0) || (nt >= ga->nr))
164 return NOTSET;
167 return ga->atomnumber[nt];
170 real get_atomtype_radius(int nt, gpp_atomtype_t ga)
172 if ((nt < 0) || (nt >= ga->nr))
174 return NOTSET;
177 return ga->radius[nt];
180 real get_atomtype_vol(int nt, gpp_atomtype_t ga)
182 if ((nt < 0) || (nt >= ga->nr))
184 return NOTSET;
187 return ga->vol[nt];
190 real get_atomtype_surftens(int nt, gpp_atomtype_t ga)
192 if ((nt < 0) || (nt >= ga->nr))
194 return NOTSET;
197 return ga->surftens[nt];
200 real get_atomtype_gb_radius(int nt, gpp_atomtype_t ga)
202 if ((nt < 0) || (nt >= ga->nr))
204 return NOTSET;
207 return ga->gb_radius[nt];
210 real get_atomtype_S_hct(int nt, gpp_atomtype_t ga)
212 if ((nt < 0) || (nt >= ga->nr))
214 return NOTSET;
217 return ga->S_hct[nt];
220 real get_atomtype_nbparam(int nt, int param, gpp_atomtype_t ga)
222 if ((nt < 0) || (nt >= ga->nr))
224 return NOTSET;
226 if ((param < 0) || (param >= MAXFORCEPARAM))
228 return NOTSET;
230 return ga->nb[nt].c[param];
233 gpp_atomtype_t init_atomtype(void)
235 gpp_atomtype_t ga;
237 snew(ga, 1);
239 ga->nr = 0;
240 ga->atom = NULL;
241 ga->atomname = NULL;
242 ga->nb = NULL;
243 ga->bondatomtype = NULL;
244 ga->radius = NULL;
245 ga->vol = NULL;
246 ga->surftens = NULL;
247 ga->atomnumber = NULL;
248 ga->gb_radius = NULL;
249 ga->S_hct = NULL;
251 return ga;
255 set_atomtype_gbparam(gpp_atomtype_t ga, int i,
256 real radius, real vol, real surftens,
257 real gb_radius, real S_hct)
259 if ( (i < 0) || (i >= ga->nr))
261 return NOTSET;
264 ga->radius[i] = radius;
265 ga->vol[i] = vol;
266 ga->surftens[i] = surftens;
267 ga->gb_radius[i] = gb_radius;
268 ga->S_hct[i] = S_hct;
270 return i;
274 int set_atomtype(int nt, gpp_atomtype_t ga, t_symtab *tab,
275 t_atom *a, const char *name, t_param *nb,
276 int bondatomtype,
277 real radius, real vol, real surftens, int atomnumber,
278 real gb_radius, real S_hct)
280 if ((nt < 0) || (nt >= ga->nr))
282 return NOTSET;
285 ga->atom[nt] = *a;
286 ga->atomname[nt] = put_symtab(tab, name);
287 ga->nb[nt] = *nb;
288 ga->bondatomtype[nt] = bondatomtype;
289 ga->radius[nt] = radius;
290 ga->vol[nt] = vol;
291 ga->surftens[nt] = surftens;
292 ga->atomnumber[nt] = atomnumber;
293 ga->gb_radius[nt] = gb_radius;
294 ga->S_hct[nt] = S_hct;
296 return nt;
299 int add_atomtype(gpp_atomtype_t ga, t_symtab *tab,
300 t_atom *a, const char *name, t_param *nb,
301 int bondatomtype,
302 real radius, real vol, real surftens, int atomnumber,
303 real gb_radius, real S_hct)
305 int i;
307 for (i = 0; (i < ga->nr); i++)
309 if (strcmp(*ga->atomname[i], name) == 0)
311 if (NULL != debug)
313 fprintf(debug, "Trying to add atomtype %s again. Skipping it.\n", name);
315 break;
318 if (i == ga->nr)
320 ga->nr++;
321 srenew(ga->atom, ga->nr);
322 srenew(ga->atomname, ga->nr);
323 srenew(ga->nb, ga->nr);
324 srenew(ga->bondatomtype, ga->nr);
325 srenew(ga->radius, ga->nr);
326 srenew(ga->vol, ga->nr);
327 srenew(ga->surftens, ga->nr);
328 srenew(ga->atomnumber, ga->nr);
329 srenew(ga->gb_radius, ga->nr);
330 srenew(ga->S_hct, ga->nr);
332 return set_atomtype(ga->nr-1, ga, tab, a, name, nb, bondatomtype, radius,
333 vol, surftens, atomnumber, gb_radius, S_hct);
335 else
337 return i;
341 void print_at (FILE * out, gpp_atomtype_t ga)
343 int i;
344 t_atom *atom = ga->atom;
345 t_param *nb = ga->nb;
347 fprintf (out, "[ %s ]\n", dir2str(d_atomtypes));
348 fprintf (out, "; %6s %8s %8s %8s %12s %12s\n",
349 "type", "mass", "charge", "particle", "c6", "c12");
350 for (i = 0; (i < ga->nr); i++)
352 fprintf(out, "%8s %8.3f %8.3f %8s %12e %12e\n",
353 *(ga->atomname[i]), atom[i].m, atom[i].q, "A",
354 nb[i].c0(), nb[i].c1());
357 fprintf (out, "\n");
360 void done_atomtype(gpp_atomtype_t ga)
362 sfree(ga->atom);
363 sfree(ga->atomname);
364 sfree(ga->nb);
365 sfree(ga->bondatomtype);
366 sfree(ga->radius);
367 sfree(ga->vol);
368 sfree(ga->gb_radius);
369 sfree(ga->S_hct);
370 sfree(ga->surftens);
371 sfree(ga->atomnumber);
372 ga->nr = 0;
373 sfree(ga);
376 static int search_atomtypes(gpp_atomtype_t ga, int *n, int typelist[],
377 int thistype,
378 t_param param[], int ftype)
380 int i, nn, nrfp, j, k, ntype, tli;
381 gmx_bool bFound = FALSE;
383 nn = *n;
384 nrfp = NRFP(ftype);
385 ntype = get_atomtype_ntypes(ga);
387 for (i = 0; (i < nn); i++)
389 if (typelist[i] == thistype)
391 /* This type number has already been added */
392 break;
395 /* Otherwise, check if the parameters are identical to any previously added type */
397 bFound = TRUE;
398 for (j = 0; j < ntype && bFound; j++)
400 /* Check nonbonded parameters */
401 for (k = 0; k < nrfp && bFound; k++)
403 bFound = (param[ntype*typelist[i]+j].c[k] == param[ntype*thistype+j].c[k]);
406 /* Check radius, volume, surftens */
407 tli = typelist[i];
408 bFound = bFound &&
409 (get_atomtype_radius(tli, ga) == get_atomtype_radius(thistype, ga)) &&
410 (get_atomtype_vol(tli, ga) == get_atomtype_vol(thistype, ga)) &&
411 (get_atomtype_surftens(tli, ga) == get_atomtype_surftens(thistype, ga)) &&
412 (get_atomtype_atomnumber(tli, ga) == get_atomtype_atomnumber(thistype, ga)) &&
413 (get_atomtype_gb_radius(tli, ga) == get_atomtype_gb_radius(thistype, ga)) &&
414 (get_atomtype_S_hct(tli, ga) == get_atomtype_S_hct(thistype, ga));
416 if (bFound)
418 break;
422 if (i == nn)
424 if (debug)
426 fprintf(debug, "Renumbering atomtype %d to %d\n", thistype, nn);
428 if (nn == ntype)
430 gmx_fatal(FARGS, "Atomtype horror n = %d, %s, %d", nn, __FILE__, __LINE__);
432 typelist[nn] = thistype;
433 nn++;
435 *n = nn;
437 return i;
440 void renum_atype(t_params plist[], gmx_mtop_t *mtop,
441 int *wall_atomtype,
442 gpp_atomtype_t ga, gmx_bool bVerbose)
444 int i, j, k, l, molt, mi, mj, nat, nrfp, ftype, ntype;
445 t_atoms *atoms;
446 t_param *nbsnew;
447 int *typelist;
448 real *new_radius;
449 real *new_vol;
450 real *new_surftens;
451 real *new_gb_radius;
452 real *new_S_hct;
453 int *new_atomnumber;
454 char ***new_atomname;
456 ntype = get_atomtype_ntypes(ga);
457 snew(typelist, ntype);
459 if (bVerbose)
461 fprintf(stderr, "renumbering atomtypes...\n");
464 /* Since the bonded interactions have been assigned now,
465 * we want to reduce the number of atom types by merging
466 * ones with identical nonbonded interactions, in addition
467 * to removing unused ones.
469 * With Generalized-Born electrostatics, or implicit solvent
470 * we also check that the atomtype radius, effective_volume
471 * and surface tension match.
473 * With QM/MM we also check that the atom numbers match
476 /* Get nonbonded interaction type */
477 if (plist[F_LJ].nr > 0)
479 ftype = F_LJ;
481 else
483 ftype = F_BHAM;
486 /* Renumber atomtypes by first making a list of which ones are actually used.
487 * We provide the list of nonbonded parameters so search_atomtypes
488 * can determine if two types should be merged.
490 nat = 0;
491 for (molt = 0; molt < mtop->nmoltype; molt++)
493 atoms = &mtop->moltype[molt].atoms;
494 for (i = 0; (i < atoms->nr); i++)
496 atoms->atom[i].type =
497 search_atomtypes(ga, &nat, typelist, atoms->atom[i].type,
498 plist[ftype].param, ftype);
499 atoms->atom[i].typeB =
500 search_atomtypes(ga, &nat, typelist, atoms->atom[i].typeB,
501 plist[ftype].param, ftype);
505 for (i = 0; i < 2; i++)
507 if (wall_atomtype[i] >= 0)
509 wall_atomtype[i] = search_atomtypes(ga, &nat, typelist, wall_atomtype[i],
510 plist[ftype].param, ftype);
514 snew(new_radius, nat);
515 snew(new_vol, nat);
516 snew(new_surftens, nat);
517 snew(new_atomnumber, nat);
518 snew(new_gb_radius, nat);
519 snew(new_S_hct, nat);
520 snew(new_atomname, nat);
522 /* We now have a list of unique atomtypes in typelist */
524 if (debug)
526 pr_ivec(debug, 0, "typelist", typelist, nat, TRUE);
529 /* Renumber nlist */
530 nbsnew = NULL;
531 snew(nbsnew, plist[ftype].nr);
533 nrfp = NRFP(ftype);
535 for (i = k = 0; (i < nat); i++)
537 mi = typelist[i];
538 for (j = 0; (j < nat); j++, k++)
540 mj = typelist[j];
541 for (l = 0; (l < nrfp); l++)
543 nbsnew[k].c[l] = plist[ftype].param[ntype*mi+mj].c[l];
546 new_radius[i] = get_atomtype_radius(mi, ga);
547 new_vol[i] = get_atomtype_vol(mi, ga);
548 new_surftens[i] = get_atomtype_surftens(mi, ga);
549 new_atomnumber[i] = get_atomtype_atomnumber(mi, ga);
550 new_gb_radius[i] = get_atomtype_gb_radius(mi, ga);
551 new_S_hct[i] = get_atomtype_S_hct(mi, ga);
552 new_atomname[i] = ga->atomname[mi];
555 for (i = 0; (i < nat*nat); i++)
557 for (l = 0; (l < nrfp); l++)
559 plist[ftype].param[i].c[l] = nbsnew[i].c[l];
562 plist[ftype].nr = i;
563 mtop->ffparams.atnr = nat;
565 sfree(ga->radius);
566 sfree(ga->vol);
567 sfree(ga->surftens);
568 sfree(ga->atomnumber);
569 sfree(ga->gb_radius);
570 sfree(ga->S_hct);
571 /* Dangling atomname pointers ? */
572 sfree(ga->atomname);
574 ga->radius = new_radius;
575 ga->vol = new_vol;
576 ga->surftens = new_surftens;
577 ga->atomnumber = new_atomnumber;
578 ga->gb_radius = new_gb_radius;
579 ga->S_hct = new_S_hct;
580 ga->atomname = new_atomname;
582 ga->nr = nat;
584 sfree(nbsnew);
585 sfree(typelist);
588 void copy_atomtype_atomtypes(gpp_atomtype_t ga, t_atomtypes *atomtypes)
590 int i, ntype;
592 /* Copy the atomtype data to the topology atomtype list */
593 ntype = get_atomtype_ntypes(ga);
594 atomtypes->nr = ntype;
595 snew(atomtypes->radius, ntype);
596 snew(atomtypes->vol, ntype);
597 snew(atomtypes->surftens, ntype);
598 snew(atomtypes->atomnumber, ntype);
599 snew(atomtypes->gb_radius, ntype);
600 snew(atomtypes->S_hct, ntype);
602 for (i = 0; i < ntype; i++)
604 atomtypes->radius[i] = ga->radius[i];
605 atomtypes->vol[i] = ga->vol[i];
606 atomtypes->surftens[i] = ga->surftens[i];
607 atomtypes->atomnumber[i] = ga->atomnumber[i];
608 atomtypes->gb_radius[i] = ga->gb_radius[i];
609 atomtypes->S_hct[i] = ga->S_hct[i];