Refactor t_pindex
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / gpp_atomtype.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) 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.
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 "gpp_atomtype.h"
40
41 #include <climits>
42 #include <cmath>
43 #include <cstring>
44
45 #include <algorithm>
46
47 #include "gromacs/gmxpreprocess/grompp_impl.h"
48 #include "gromacs/gmxpreprocess/notset.h"
49 #include "gromacs/gmxpreprocess/topdirs.h"
50 #include "gromacs/gmxpreprocess/toputil.h"
51 #include "gromacs/math/vecdump.h"
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/topology/symtab.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/smalloc.h"
57
58 struct AtomTypeData
59 {
60     //! Explicit constructor.
61     AtomTypeData(const t_atom            &a,
62                  char                   **name,
63                  const InteractionOfType &nb,
64                  const int                bondAtomType,
65                  const int                atomNumber) :
66         atom_(a), name_(name), nb_(nb),
67         bondAtomType_(bondAtomType),
68         atomNumber_(atomNumber)
69     {
70     }
71     //! Actual atom data.
72     t_atom             atom_;
73     //! Atom name.
74     char             **name_;
75     //! Nonbonded data.
76     InteractionOfType  nb_;
77     //! Bonded atomtype for the type.
78     int                bondAtomType_;
79     //! Atom number for the atom type.
80     int                atomNumber_;
81 };
82
83 class PreprocessingAtomTypes::Impl
84 {
85     public:
86         //! The number for currently loaded entries.
87         size_t size() const { return types.size(); }
88         //! The actual atom type data.
89         std::vector<AtomTypeData> types;
90 };
91
92 bool PreprocessingAtomTypes::isSet(int nt) const
93 {
94     return ((nt >= 0) && (nt < gmx::ssize(*this)));
95 }
96
97 int PreprocessingAtomTypes::atomTypeFromName(const std::string &str) const
98 {
99     /* Atom types are always case sensitive */
100     auto found = std::find_if(impl_->types.begin(), impl_->types.end(),
101                               [&str](const auto &type)
102                               { return str == *type.name_; });
103     if (found == impl_->types.end())
104     {
105         return NOTSET;
106     }
107     else
108     {
109         return std::distance(impl_->types.begin(), found);
110     }
111 }
112
113 size_t PreprocessingAtomTypes::size() const
114 {
115     return impl_->size();
116 }
117
118 const char *PreprocessingAtomTypes::atomNameFromAtomType(int nt) const
119 {
120     return isSet(nt) ? *(impl_->types[nt].name_) : nullptr;
121 }
122
123 real PreprocessingAtomTypes::atomMassFromAtomType(int nt) const
124 {
125     return isSet(nt) ? impl_->types[nt].atom_.m : NOTSET;
126 }
127
128 real PreprocessingAtomTypes::atomChargeFromAtomType(int nt) const
129 {
130     return isSet(nt) ? impl_->types[nt].atom_.q : NOTSET;
131 }
132
133 int PreprocessingAtomTypes::atomParticleTypeFromAtomType(int nt) const
134 {
135     return isSet(nt) ? impl_->types[nt].atom_.ptype : NOTSET;
136 }
137
138 int PreprocessingAtomTypes::bondAtomTypeFromAtomType(int nt) const
139 {
140     return isSet(nt) ? impl_->types[nt].bondAtomType_ : NOTSET;
141 }
142
143 int PreprocessingAtomTypes::atomNumberFromAtomType(int nt) const
144 {
145     return isSet(nt) ? impl_->types[nt].atomNumber_ : NOTSET;
146 }
147
148 real PreprocessingAtomTypes::atomNonBondedParamFromAtomType(int nt, int param) const
149 {
150     if (!isSet(nt))
151     {
152         return NOTSET;
153     }
154     gmx::ArrayRef<const real> forceParam = impl_->types[nt].nb_.forceParam();
155     if ((param < 0) || (param >= MAXFORCEPARAM))
156     {
157         return NOTSET;
158     }
159     return forceParam[param];
160 }
161
162 PreprocessingAtomTypes::PreprocessingAtomTypes()
163     : impl_(new Impl)
164 {}
165
166 PreprocessingAtomTypes::PreprocessingAtomTypes(PreprocessingAtomTypes &&old) noexcept
167     : impl_(std::move(old.impl_))
168 {}
169
170 PreprocessingAtomTypes &PreprocessingAtomTypes::operator=(PreprocessingAtomTypes &&old) noexcept
171 {
172     impl_ = std::move(old.impl_);
173     return *this;
174 }
175
176 PreprocessingAtomTypes::~PreprocessingAtomTypes()
177 {}
178
179 int PreprocessingAtomTypes::addType(t_symtab                *tab,
180                                     const t_atom            &a,
181                                     const std::string       &name,
182                                     const InteractionOfType &nb,
183                                     int                      bondAtomType,
184                                     int                      atomNumber)
185 {
186     int position = atomTypeFromName(name);
187     if (position == NOTSET)
188     {
189         impl_->types.emplace_back(a,
190                                   put_symtab(tab, name.c_str()),
191                                   nb,
192                                   bondAtomType,
193                                   atomNumber);
194         return atomTypeFromName(name);
195     }
196     else
197     {
198         return position;
199     }
200 }
201
202 int PreprocessingAtomTypes::setType(int                      nt,
203                                     t_symtab                *tab,
204                                     const t_atom            &a,
205                                     const std::string       &name,
206                                     const InteractionOfType &nb,
207                                     int                      bondAtomType,
208                                     int                      atomNumber)
209 {
210     if (!isSet(nt))
211     {
212         return NOTSET;
213     }
214
215     impl_->types[nt].atom_         = a;
216     impl_->types[nt].name_         = put_symtab(tab, name.c_str());
217     impl_->types[nt].nb_           = nb;
218     impl_->types[nt].bondAtomType_ = bondAtomType;
219     impl_->types[nt].atomNumber_   = atomNumber;
220
221     return nt;
222 }
223
224 void PreprocessingAtomTypes::printTypes(FILE * out)
225 {
226     fprintf (out, "[ %s ]\n", dir2str(Directive::d_atomtypes));
227     fprintf (out, "; %6s  %8s  %8s  %8s  %12s  %12s\n",
228              "type", "mass", "charge", "particle", "c6", "c12");
229     for (auto &entry : impl_->types)
230     {
231         fprintf(out, "%8s  %8.3f  %8.3f  %8s  %12e  %12e\n",
232                 *(entry.name_), entry.atom_.m, entry.atom_.q, "A",
233                 entry.nb_.c0(), entry.nb_.c1());
234     }
235
236     fprintf (out, "\n");
237 }
238
239 static int search_atomtypes(const PreprocessingAtomTypes          *ga,
240                             int                                   *n,
241                             gmx::ArrayRef<int>                     typelist,
242                             int                                    thistype,
243                             gmx::ArrayRef<const InteractionOfType> interactionTypes,
244                             int                                    ftype)
245 {
246     int nn    = *n;
247     int nrfp  = NRFP(ftype);
248     int ntype = ga->size();
249
250     int i;
251     for (i = 0; (i < nn); i++)
252     {
253         if (typelist[i] == thistype)
254         {
255             /* This type number has already been added */
256             break;
257         }
258
259         /* Otherwise, check if the parameters are identical to any previously added type */
260
261         bool bFound = true;
262         for (int j = 0; j < ntype && bFound; j++)
263         {
264             /* Check nonbonded parameters */
265             gmx::ArrayRef<const real> forceParam1 = interactionTypes[ntype*typelist[i]+j].forceParam();
266             gmx::ArrayRef<const real> forceParam2 = interactionTypes[ntype*thistype+j].forceParam();
267             for (int k = 0; (k < nrfp) && bFound; k++)
268             {
269                 bFound = forceParam1[k] == forceParam2[k];
270             }
271
272             /* Check atomnumber */
273             int tli    = typelist[i];
274             bFound = bFound &&
275                 (ga->atomNumberFromAtomType(tli) == ga->atomNumberFromAtomType(thistype));
276         }
277         if (bFound)
278         {
279             break;
280         }
281     }
282
283     if (i == nn)
284     {
285         if (nn == ntype)
286         {
287             gmx_fatal(FARGS, "Atomtype horror n = %d, %s, %d", nn, __FILE__, __LINE__);
288         }
289         typelist[nn] = thistype;
290         nn++;
291     }
292     *n = nn;
293
294     return i;
295 }
296
297 void PreprocessingAtomTypes::renumberTypes(gmx::ArrayRef<InteractionsOfType>        plist,
298                                            gmx_mtop_t                              *mtop,
299                                            int                                     *wall_atomtype,
300                                            bool                                     bVerbose)
301 {
302     int         nat, ftype, ntype;
303
304     ntype = size();
305     std::vector<int> typelist(ntype);
306
307     if (bVerbose)
308     {
309         fprintf(stderr, "renumbering atomtypes...\n");
310     }
311
312     /* Since the bonded interactions have been assigned now,
313      * we want to reduce the number of atom types by merging
314      * ones with identical nonbonded interactions, in addition
315      * to removing unused ones.
316      *
317      * With QM/MM we also check that the atom numbers match
318      */
319
320     /* Get nonbonded interaction type */
321     if (plist[F_LJ].size() > 0)
322     {
323         ftype = F_LJ;
324     }
325     else
326     {
327         ftype = F_BHAM;
328     }
329
330     /* Renumber atomtypes by first making a list of which ones are actually used.
331      * We provide the list of nonbonded parameters so search_atomtypes
332      * can determine if two types should be merged.
333      */
334     nat = 0;
335     for (const gmx_moltype_t &moltype : mtop->moltype)
336     {
337         const t_atoms *atoms = &moltype.atoms;
338         for (int i = 0; (i < atoms->nr); i++)
339         {
340             atoms->atom[i].type =
341                 search_atomtypes(this, &nat, typelist, atoms->atom[i].type,
342                                  plist[ftype].interactionTypes, ftype);
343             atoms->atom[i].typeB =
344                 search_atomtypes(this, &nat, typelist, atoms->atom[i].typeB,
345                                  plist[ftype].interactionTypes, ftype);
346         }
347     }
348
349     for (int i = 0; i < 2; i++)
350     {
351         if (wall_atomtype[i] >= 0)
352         {
353             wall_atomtype[i] = search_atomtypes(this, &nat, typelist, wall_atomtype[i],
354                                                 plist[ftype].interactionTypes, ftype);
355         }
356     }
357
358     std::vector<AtomTypeData> new_types;
359     /* We now have a list of unique atomtypes in typelist */
360
361     /* Renumber nlist */
362     std::vector<InteractionOfType> nbsnew;
363
364     for (int i = 0; (i < nat); i++)
365     {
366         int mi = typelist[i];
367         for (int j = 0; (j < nat); j++)
368         {
369             int                      mj              = typelist[j];
370             const InteractionOfType &interactionType = plist[ftype].interactionTypes[ntype*mi+mj];
371             nbsnew.emplace_back(interactionType.atoms(), interactionType.forceParam(), interactionType.interactionTypeName());
372         }
373         new_types.push_back(impl_->types[mi]);
374     }
375
376     mtop->ffparams.atnr = nat;
377
378     impl_->types                  = new_types;
379     plist[ftype].interactionTypes = nbsnew;
380 }
381
382 void PreprocessingAtomTypes::copyTot_atomtypes(t_atomtypes *atomtypes) const
383 {
384     /* Copy the atomtype data to the topology atomtype list */
385     int ntype         = size();
386     atomtypes->nr = ntype;
387     snew(atomtypes->atomnumber, ntype);
388
389     for (int i = 0; i < ntype; i++)
390     {
391         atomtypes->atomnumber[i] = impl_->types[i].atomNumber_;
392     }
393 }