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,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.
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.
39 #include "gpp_atomtype.h"
46 #include "gromacs/gmxpreprocess/notset.h"
47 #include "gromacs/gmxpreprocess/topdirs.h"
48 #include "gromacs/gmxpreprocess/toputil.h"
49 #include "gromacs/math/vecdump.h"
50 #include "gromacs/topology/ifunc.h"
51 #include "gromacs/topology/symtab.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/smalloc.h"
56 typedef struct gpp_atomtype {
57 int nr; /* The number of atomtypes */
58 t_atom *atom; /* Array of atoms */
59 char ***atomname; /* Names of the atomtypes */
60 t_param *nb; /* Nonbonded force default params */
61 int *bondatomtype; /* The bond_atomtype for each atomtype */
62 int *atomnumber; /* Atomic number, used for QM/MM */
65 int get_atomtype_type(const char *str, gpp_atomtype_t ga)
69 /* Atom types are always case sensitive */
70 for (i = 0; (i < ga->nr); i++)
72 if (strcmp(str, *(ga->atomname[i])) == 0)
81 int get_atomtype_ntypes(gpp_atomtype_t ga)
86 char *get_atomtype_name(int nt, gpp_atomtype_t ga)
88 if ((nt < 0) || (nt >= ga->nr))
93 return *(ga->atomname[nt]);
96 real get_atomtype_massA(int nt, gpp_atomtype_t ga)
98 if ((nt < 0) || (nt >= ga->nr))
103 return ga->atom[nt].m;
106 real get_atomtype_massB(int nt, gpp_atomtype_t ga)
108 if ((nt < 0) || (nt >= ga->nr))
113 return ga->atom[nt].mB;
116 real get_atomtype_qA(int nt, gpp_atomtype_t ga)
118 if ((nt < 0) || (nt >= ga->nr))
123 return ga->atom[nt].q;
126 real get_atomtype_qB(int nt, gpp_atomtype_t ga)
128 if ((nt < 0) || (nt >= ga->nr))
133 return ga->atom[nt].qB;
136 int get_atomtype_ptype(int nt, gpp_atomtype_t ga)
138 if ((nt < 0) || (nt >= ga->nr))
143 return ga->atom[nt].ptype;
146 int get_atomtype_batype(int nt, const gpp_atomtype* ga)
148 if ((nt < 0) || (nt >= ga->nr))
153 return ga->bondatomtype[nt];
156 int get_atomtype_atomnumber(int nt, gpp_atomtype_t ga)
158 if ((nt < 0) || (nt >= ga->nr))
163 return ga->atomnumber[nt];
166 real get_atomtype_nbparam(int nt, int param, gpp_atomtype_t ga)
168 if ((nt < 0) || (nt >= ga->nr))
172 if ((param < 0) || (param >= MAXFORCEPARAM))
176 return ga->nb[nt].c[param];
179 gpp_atomtype_t init_atomtype(void)
187 ga->atomname = nullptr;
189 ga->bondatomtype = nullptr;
190 ga->atomnumber = nullptr;
195 int set_atomtype(int nt, gpp_atomtype_t ga, t_symtab *tab,
196 t_atom *a, const char *name, t_param *nb,
197 int bondatomtype, int atomnumber)
199 if ((nt < 0) || (nt >= ga->nr))
205 ga->atomname[nt] = put_symtab(tab, name);
207 ga->bondatomtype[nt] = bondatomtype;
208 ga->atomnumber[nt] = atomnumber;
213 int add_atomtype(gpp_atomtype_t ga, t_symtab *tab,
214 t_atom *a, const char *name, t_param *nb,
215 int bondatomtype, int atomnumber)
219 for (i = 0; (i < ga->nr); i++)
221 if (strcmp(*ga->atomname[i], name) == 0)
229 srenew(ga->atom, ga->nr);
230 srenew(ga->atomname, ga->nr);
231 srenew(ga->nb, ga->nr);
232 srenew(ga->bondatomtype, ga->nr);
233 srenew(ga->atomnumber, ga->nr);
235 return set_atomtype(ga->nr-1, ga, tab, a, name, nb, bondatomtype, atomnumber);
243 void print_at (FILE * out, gpp_atomtype_t ga)
246 t_atom *atom = ga->atom;
247 t_param *nb = ga->nb;
249 fprintf (out, "[ %s ]\n", dir2str(d_atomtypes));
250 fprintf (out, "; %6s %8s %8s %8s %12s %12s\n",
251 "type", "mass", "charge", "particle", "c6", "c12");
252 for (i = 0; (i < ga->nr); i++)
254 fprintf(out, "%8s %8.3f %8.3f %8s %12e %12e\n",
255 *(ga->atomname[i]), atom[i].m, atom[i].q, "A",
256 nb[i].c0(), nb[i].c1());
262 void done_atomtype(gpp_atomtype_t ga)
267 sfree(ga->bondatomtype);
268 sfree(ga->atomnumber);
273 static int search_atomtypes(gpp_atomtype_t ga, int *n, int typelist[],
275 t_param param[], int ftype)
277 int i, nn, nrfp, j, k, ntype, tli;
282 ntype = get_atomtype_ntypes(ga);
284 for (i = 0; (i < nn); i++)
286 if (typelist[i] == thistype)
288 /* This type number has already been added */
292 /* Otherwise, check if the parameters are identical to any previously added type */
295 for (j = 0; j < ntype && bFound; j++)
297 /* Check nonbonded parameters */
298 for (k = 0; k < nrfp && bFound; k++)
300 bFound = (param[ntype*typelist[i]+j].c[k] == param[ntype*thistype+j].c[k]);
303 /* Check atomnumber */
306 (get_atomtype_atomnumber(tli, ga) == get_atomtype_atomnumber(thistype, ga));
318 gmx_fatal(FARGS, "Atomtype horror n = %d, %s, %d", nn, __FILE__, __LINE__);
320 typelist[nn] = thistype;
328 void renum_atype(t_params plist[], gmx_mtop_t *mtop,
330 gpp_atomtype_t ga, bool bVerbose)
332 int i, j, k, l, mi, mj, nat, nrfp, ftype, ntype;
337 char ***new_atomname;
339 ntype = get_atomtype_ntypes(ga);
340 snew(typelist, ntype);
344 fprintf(stderr, "renumbering atomtypes...\n");
347 /* Since the bonded interactions have been assigned now,
348 * we want to reduce the number of atom types by merging
349 * ones with identical nonbonded interactions, in addition
350 * to removing unused ones.
352 * With QM/MM we also check that the atom numbers match
355 /* Get nonbonded interaction type */
356 if (plist[F_LJ].nr > 0)
365 /* Renumber atomtypes by first making a list of which ones are actually used.
366 * We provide the list of nonbonded parameters so search_atomtypes
367 * can determine if two types should be merged.
370 for (gmx_moltype_t &moltype : mtop->moltype)
372 atoms = &moltype.atoms;
373 for (i = 0; (i < atoms->nr); i++)
375 atoms->atom[i].type =
376 search_atomtypes(ga, &nat, typelist, atoms->atom[i].type,
377 plist[ftype].param, ftype);
378 atoms->atom[i].typeB =
379 search_atomtypes(ga, &nat, typelist, atoms->atom[i].typeB,
380 plist[ftype].param, ftype);
384 for (i = 0; i < 2; i++)
386 if (wall_atomtype[i] >= 0)
388 wall_atomtype[i] = search_atomtypes(ga, &nat, typelist, wall_atomtype[i],
389 plist[ftype].param, ftype);
393 snew(new_atomnumber, nat);
394 snew(new_atomname, nat);
395 /* We now have a list of unique atomtypes in typelist */
399 snew(nbsnew, plist[ftype].nr);
403 for (i = k = 0; (i < nat); i++)
406 for (j = 0; (j < nat); j++, k++)
409 for (l = 0; (l < nrfp); l++)
411 nbsnew[k].c[l] = plist[ftype].param[ntype*mi+mj].c[l];
414 new_atomnumber[i] = get_atomtype_atomnumber(mi, ga);
415 new_atomname[i] = ga->atomname[mi];
418 for (i = 0; (i < nat*nat); i++)
420 for (l = 0; (l < nrfp); l++)
422 plist[ftype].param[i].c[l] = nbsnew[i].c[l];
426 mtop->ffparams.atnr = nat;
428 sfree(ga->atomnumber);
429 /* Dangling atomname pointers ? */
432 ga->atomnumber = new_atomnumber;
433 ga->atomname = new_atomname;
441 void copy_atomtype_atomtypes(gpp_atomtype_t ga, t_atomtypes *atomtypes)
445 /* Copy the atomtype data to the topology atomtype list */
446 ntype = get_atomtype_ntypes(ga);
447 atomtypes->nr = ntype;
448 snew(atomtypes->atomnumber, ntype);
450 for (i = 0; i < ntype; i++)
452 atomtypes->atomnumber[i] = ga->atomnumber[i];