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,2019, 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"
45 #include "gromacs/gmxpreprocess/notset.h"
46 #include "gromacs/gmxpreprocess/topdirs.h"
47 #include "gromacs/gmxpreprocess/toputil.h"
48 #include "gromacs/math/vecdump.h"
49 #include "gromacs/topology/ifunc.h"
50 #include "gromacs/topology/symtab.h"
51 #include "gromacs/topology/topology.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 int *atomnumber; /* Atomic number, used for QM/MM */
64 int get_atomtype_type(const char *str, gpp_atomtype_t ga)
68 /* Atom types are always case sensitive */
69 for (i = 0; (i < ga->nr); i++)
71 if (strcmp(str, *(ga->atomname[i])) == 0)
80 int get_atomtype_ntypes(gpp_atomtype_t ga)
85 char *get_atomtype_name(int nt, gpp_atomtype_t ga)
87 if ((nt < 0) || (nt >= ga->nr))
92 return *(ga->atomname[nt]);
95 real get_atomtype_massA(int nt, gpp_atomtype_t ga)
97 if ((nt < 0) || (nt >= ga->nr))
102 return ga->atom[nt].m;
105 real get_atomtype_massB(int nt, gpp_atomtype_t ga)
107 if ((nt < 0) || (nt >= ga->nr))
112 return ga->atom[nt].mB;
115 real get_atomtype_qA(int nt, gpp_atomtype_t ga)
117 if ((nt < 0) || (nt >= ga->nr))
122 return ga->atom[nt].q;
125 real get_atomtype_qB(int nt, gpp_atomtype_t ga)
127 if ((nt < 0) || (nt >= ga->nr))
132 return ga->atom[nt].qB;
135 int get_atomtype_ptype(int nt, gpp_atomtype_t ga)
137 if ((nt < 0) || (nt >= ga->nr))
142 return ga->atom[nt].ptype;
145 int get_atomtype_batype(int nt, const gpp_atomtype* ga)
147 if ((nt < 0) || (nt >= ga->nr))
152 return ga->bondatomtype[nt];
155 int get_atomtype_atomnumber(int nt, gpp_atomtype_t ga)
157 if ((nt < 0) || (nt >= ga->nr))
162 return ga->atomnumber[nt];
165 real get_atomtype_nbparam(int nt, int param, gpp_atomtype_t ga)
167 if ((nt < 0) || (nt >= ga->nr))
171 if ((param < 0) || (param >= MAXFORCEPARAM))
175 return ga->nb[nt].c[param];
178 gpp_atomtype_t init_atomtype()
186 ga->atomname = nullptr;
188 ga->bondatomtype = nullptr;
189 ga->atomnumber = nullptr;
194 int set_atomtype(int nt, gpp_atomtype_t ga, t_symtab *tab,
195 t_atom *a, const char *name, t_param *nb,
196 int bondatomtype, int atomnumber)
198 if ((nt < 0) || (nt >= ga->nr))
204 ga->atomname[nt] = put_symtab(tab, name);
206 ga->bondatomtype[nt] = bondatomtype;
207 ga->atomnumber[nt] = atomnumber;
212 int add_atomtype(gpp_atomtype_t ga, t_symtab *tab,
213 t_atom *a, const char *name, t_param *nb,
214 int bondatomtype, int atomnumber)
218 for (i = 0; (i < ga->nr); i++)
220 if (strcmp(*ga->atomname[i], name) == 0)
228 srenew(ga->atom, ga->nr);
229 srenew(ga->atomname, ga->nr);
230 srenew(ga->nb, ga->nr);
231 srenew(ga->bondatomtype, ga->nr);
232 srenew(ga->atomnumber, ga->nr);
234 return set_atomtype(ga->nr-1, ga, tab, a, name, nb, bondatomtype, atomnumber);
242 void print_at (FILE * out, gpp_atomtype_t ga)
245 t_atom *atom = ga->atom;
246 t_param *nb = ga->nb;
248 fprintf (out, "[ %s ]\n", dir2str(Directive::d_atomtypes));
249 fprintf (out, "; %6s %8s %8s %8s %12s %12s\n",
250 "type", "mass", "charge", "particle", "c6", "c12");
251 for (i = 0; (i < ga->nr); i++)
253 fprintf(out, "%8s %8.3f %8.3f %8s %12e %12e\n",
254 *(ga->atomname[i]), atom[i].m, atom[i].q, "A",
255 nb[i].c0(), nb[i].c1());
261 void done_atomtype(gpp_atomtype_t ga)
266 sfree(ga->bondatomtype);
267 sfree(ga->atomnumber);
272 static int search_atomtypes(gpp_atomtype_t ga, int *n, int typelist[],
274 t_param param[], int ftype)
276 int i, nn, nrfp, j, k, ntype, tli;
281 ntype = get_atomtype_ntypes(ga);
283 for (i = 0; (i < nn); i++)
285 if (typelist[i] == thistype)
287 /* This type number has already been added */
291 /* Otherwise, check if the parameters are identical to any previously added type */
294 for (j = 0; j < ntype && bFound; j++)
296 /* Check nonbonded parameters */
297 for (k = 0; k < nrfp && bFound; k++)
299 bFound = (param[ntype*typelist[i]+j].c[k] == param[ntype*thistype+j].c[k]);
302 /* Check atomnumber */
305 (get_atomtype_atomnumber(tli, ga) == get_atomtype_atomnumber(thistype, ga));
317 gmx_fatal(FARGS, "Atomtype horror n = %d, %s, %d", nn, __FILE__, __LINE__);
319 typelist[nn] = thistype;
327 void renum_atype(t_params plist[], gmx_mtop_t *mtop,
329 gpp_atomtype_t ga, bool bVerbose)
331 int i, j, k, l, mi, mj, nat, nrfp, ftype, ntype;
336 char ***new_atomname;
338 ntype = get_atomtype_ntypes(ga);
339 snew(typelist, ntype);
343 fprintf(stderr, "renumbering atomtypes...\n");
346 /* Since the bonded interactions have been assigned now,
347 * we want to reduce the number of atom types by merging
348 * ones with identical nonbonded interactions, in addition
349 * to removing unused ones.
351 * With QM/MM we also check that the atom numbers match
354 /* Get nonbonded interaction type */
355 if (plist[F_LJ].nr > 0)
364 /* Renumber atomtypes by first making a list of which ones are actually used.
365 * We provide the list of nonbonded parameters so search_atomtypes
366 * can determine if two types should be merged.
369 for (gmx_moltype_t &moltype : mtop->moltype)
371 atoms = &moltype.atoms;
372 for (i = 0; (i < atoms->nr); i++)
374 atoms->atom[i].type =
375 search_atomtypes(ga, &nat, typelist, atoms->atom[i].type,
376 plist[ftype].param, ftype);
377 atoms->atom[i].typeB =
378 search_atomtypes(ga, &nat, typelist, atoms->atom[i].typeB,
379 plist[ftype].param, ftype);
383 for (i = 0; i < 2; i++)
385 if (wall_atomtype[i] >= 0)
387 wall_atomtype[i] = search_atomtypes(ga, &nat, typelist, wall_atomtype[i],
388 plist[ftype].param, ftype);
392 snew(new_atomnumber, nat);
393 snew(new_atomname, nat);
394 /* We now have a list of unique atomtypes in typelist */
398 snew(nbsnew, plist[ftype].nr);
402 for (i = k = 0; (i < nat); i++)
405 for (j = 0; (j < nat); j++, k++)
408 for (l = 0; (l < nrfp); l++)
410 nbsnew[k].c[l] = plist[ftype].param[ntype*mi+mj].c[l];
413 new_atomnumber[i] = get_atomtype_atomnumber(mi, ga);
414 new_atomname[i] = ga->atomname[mi];
417 for (i = 0; (i < nat*nat); i++)
419 for (l = 0; (l < nrfp); l++)
421 plist[ftype].param[i].c[l] = nbsnew[i].c[l];
425 mtop->ffparams.atnr = nat;
427 sfree(ga->atomnumber);
428 /* Dangling atomname pointers ? */
431 ga->atomnumber = new_atomnumber;
432 ga->atomname = new_atomname;
440 void copy_atomtype_atomtypes(gpp_atomtype_t ga, t_atomtypes *atomtypes)
444 /* Copy the atomtype data to the topology atomtype list */
445 ntype = get_atomtype_ntypes(ga);
446 atomtypes->nr = ntype;
447 snew(atomtypes->atomnumber, ntype);
449 for (i = 0; i < ntype; i++)
451 atomtypes->atomnumber[i] = ga->atomnumber[i];