Apply re-formatting to C++ in src/ tree.
[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 by the GROMACS development team.
7  * Copyright (c) 2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "gpp_atomtype.h"
41
42 #include <climits>
43 #include <cmath>
44 #include <cstring>
45
46 #include <algorithm>
47
48 #include "gromacs/gmxpreprocess/grompp_impl.h"
49 #include "gromacs/gmxpreprocess/notset.h"
50 #include "gromacs/gmxpreprocess/topdirs.h"
51 #include "gromacs/gmxpreprocess/toputil.h"
52 #include "gromacs/math/vecdump.h"
53 #include "gromacs/topology/ifunc.h"
54 #include "gromacs/topology/symtab.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/smalloc.h"
58
59 struct AtomTypeData
60 {
61     //! Explicit constructor.
62     AtomTypeData(const t_atom& a, char** name, const InteractionOfType& nb, const int bondAtomType, const int atomNumber) :
63         atom_(a),
64         name_(name),
65         nb_(nb),
66         bondAtomType_(bondAtomType),
67         atomNumber_(atomNumber)
68     {
69     }
70     //! Actual atom data.
71     t_atom atom_;
72     //! Atom name.
73     char** name_;
74     //! Nonbonded data.
75     InteractionOfType nb_;
76     //! Bonded atomtype for the type.
77     int bondAtomType_;
78     //! Atom number for the atom type.
79     int atomNumber_;
80 };
81
82 class PreprocessingAtomTypes::Impl
83 {
84 public:
85     //! The number for currently loaded entries.
86     size_t size() const { return types.size(); }
87     //! The actual atom type data.
88     std::vector<AtomTypeData> types;
89 };
90
91 bool PreprocessingAtomTypes::isSet(int nt) const
92 {
93     return ((nt >= 0) && (nt < gmx::ssize(*this)));
94 }
95
96 int PreprocessingAtomTypes::atomTypeFromName(const std::string& str) const
97 {
98     /* Atom types are always case sensitive */
99     auto found = std::find_if(impl_->types.begin(), impl_->types.end(), [&str](const auto& type) {
100         return str == *type.name_;
101     });
102     if (found == impl_->types.end())
103     {
104         return NOTSET;
105     }
106     else
107     {
108         return std::distance(impl_->types.begin(), found);
109     }
110 }
111
112 size_t PreprocessingAtomTypes::size() const
113 {
114     return impl_->size();
115 }
116
117 const char* PreprocessingAtomTypes::atomNameFromAtomType(int nt) const
118 {
119     return isSet(nt) ? *(impl_->types[nt].name_) : nullptr;
120 }
121
122 real PreprocessingAtomTypes::atomMassFromAtomType(int nt) const
123 {
124     return isSet(nt) ? impl_->types[nt].atom_.m : NOTSET;
125 }
126
127 real PreprocessingAtomTypes::atomChargeFromAtomType(int nt) const
128 {
129     return isSet(nt) ? impl_->types[nt].atom_.q : NOTSET;
130 }
131
132 int PreprocessingAtomTypes::atomParticleTypeFromAtomType(int nt) const
133 {
134     return isSet(nt) ? impl_->types[nt].atom_.ptype : NOTSET;
135 }
136
137 int PreprocessingAtomTypes::bondAtomTypeFromAtomType(int nt) const
138 {
139     return isSet(nt) ? impl_->types[nt].bondAtomType_ : NOTSET;
140 }
141
142 int PreprocessingAtomTypes::atomNumberFromAtomType(int nt) const
143 {
144     return isSet(nt) ? impl_->types[nt].atomNumber_ : NOTSET;
145 }
146
147 real PreprocessingAtomTypes::atomNonBondedParamFromAtomType(int nt, int param) const
148 {
149     if (!isSet(nt))
150     {
151         return NOTSET;
152     }
153     gmx::ArrayRef<const real> forceParam = impl_->types[nt].nb_.forceParam();
154     if ((param < 0) || (param >= MAXFORCEPARAM))
155     {
156         return NOTSET;
157     }
158     return forceParam[param];
159 }
160
161 PreprocessingAtomTypes::PreprocessingAtomTypes() : impl_(new Impl) {}
162
163 PreprocessingAtomTypes::PreprocessingAtomTypes(PreprocessingAtomTypes&& old) noexcept :
164     impl_(std::move(old.impl_))
165 {
166 }
167
168 PreprocessingAtomTypes& PreprocessingAtomTypes::operator=(PreprocessingAtomTypes&& old) noexcept
169 {
170     impl_ = std::move(old.impl_);
171     return *this;
172 }
173
174 PreprocessingAtomTypes::~PreprocessingAtomTypes() {}
175
176 int PreprocessingAtomTypes::addType(t_symtab*                tab,
177                                     const t_atom&            a,
178                                     const std::string&       name,
179                                     const InteractionOfType& nb,
180                                     int                      bondAtomType,
181                                     int                      atomNumber)
182 {
183     int position = atomTypeFromName(name);
184     if (position == NOTSET)
185     {
186         impl_->types.emplace_back(a, put_symtab(tab, name.c_str()), nb, bondAtomType, atomNumber);
187         return atomTypeFromName(name);
188     }
189     else
190     {
191         return position;
192     }
193 }
194
195 int PreprocessingAtomTypes::setType(int                      nt,
196                                     t_symtab*                tab,
197                                     const t_atom&            a,
198                                     const std::string&       name,
199                                     const InteractionOfType& nb,
200                                     int                      bondAtomType,
201                                     int                      atomNumber)
202 {
203     if (!isSet(nt))
204     {
205         return NOTSET;
206     }
207
208     impl_->types[nt].atom_         = a;
209     impl_->types[nt].name_         = put_symtab(tab, name.c_str());
210     impl_->types[nt].nb_           = nb;
211     impl_->types[nt].bondAtomType_ = bondAtomType;
212     impl_->types[nt].atomNumber_   = atomNumber;
213
214     return nt;
215 }
216
217 void PreprocessingAtomTypes::printTypes(FILE* out)
218 {
219     fprintf(out, "[ %s ]\n", dir2str(Directive::d_atomtypes));
220     fprintf(out,
221             "; %6s  %8s  %8s  %8s  %12s  %12s\n",
222             "type",
223             "mass",
224             "charge",
225             "particle",
226             "c6",
227             "c12");
228     for (auto& entry : impl_->types)
229     {
230         fprintf(out,
231                 "%8s  %8.3f  %8.3f  %8s  %12e  %12e\n",
232                 *(entry.name_),
233                 entry.atom_.m,
234                 entry.atom_.q,
235                 "A",
236                 entry.nb_.c0(),
237                 entry.nb_.c1());
238     }
239
240     fprintf(out, "\n");
241 }
242
243 static int search_atomtypes(const PreprocessingAtomTypes*          ga,
244                             int*                                   n,
245                             gmx::ArrayRef<int>                     typelist,
246                             int                                    thistype,
247                             gmx::ArrayRef<const InteractionOfType> interactionTypes,
248                             int                                    ftype)
249 {
250     int nn    = *n;
251     int nrfp  = NRFP(ftype);
252     int ntype = ga->size();
253
254     int i;
255     for (i = 0; (i < nn); i++)
256     {
257         if (typelist[i] == thistype)
258         {
259             /* This type number has already been added */
260             break;
261         }
262
263         /* Otherwise, check if the parameters are identical to any previously added type */
264
265         bool bFound = true;
266         for (int j = 0; j < ntype && bFound; j++)
267         {
268             /* Check nonbonded parameters */
269             gmx::ArrayRef<const real> forceParam1 =
270                     interactionTypes[ntype * typelist[i] + j].forceParam();
271             gmx::ArrayRef<const real> forceParam2 = interactionTypes[ntype * thistype + j].forceParam();
272             for (int k = 0; (k < nrfp) && bFound; k++)
273             {
274                 bFound = forceParam1[k] == forceParam2[k];
275             }
276
277             /* Check atomnumber */
278             int tli = typelist[i];
279             bFound = bFound && (ga->atomNumberFromAtomType(tli) == ga->atomNumberFromAtomType(thistype));
280         }
281         if (bFound)
282         {
283             break;
284         }
285     }
286
287     if (i == nn)
288     {
289         if (nn == ntype)
290         {
291             gmx_fatal(FARGS, "Atomtype horror n = %d, %s, %d", nn, __FILE__, __LINE__);
292         }
293         typelist[nn] = thistype;
294         nn++;
295     }
296     *n = nn;
297
298     return i;
299 }
300
301 void PreprocessingAtomTypes::renumberTypes(gmx::ArrayRef<InteractionsOfType> plist,
302                                            gmx_mtop_t*                       mtop,
303                                            int*                              wall_atomtype,
304                                            bool                              bVerbose)
305 {
306     int nat, ftype, ntype;
307
308     ntype = size();
309     std::vector<int> typelist(ntype);
310
311     if (bVerbose)
312     {
313         fprintf(stderr, "renumbering atomtypes...\n");
314     }
315
316     /* Since the bonded interactions have been assigned now,
317      * we want to reduce the number of atom types by merging
318      * ones with identical nonbonded interactions, in addition
319      * to removing unused ones.
320      *
321      * With QM/MM we also check that the atom numbers match
322      */
323
324     /* Get nonbonded interaction type */
325     if (plist[F_LJ].size() > 0)
326     {
327         ftype = F_LJ;
328     }
329     else
330     {
331         ftype = F_BHAM;
332     }
333
334     /* Renumber atomtypes by first making a list of which ones are actually used.
335      * We provide the list of nonbonded parameters so search_atomtypes
336      * can determine if two types should be merged.
337      */
338     nat = 0;
339     for (const gmx_moltype_t& moltype : mtop->moltype)
340     {
341         const t_atoms* atoms = &moltype.atoms;
342         for (int i = 0; (i < atoms->nr); i++)
343         {
344             atoms->atom[i].type = search_atomtypes(
345                     this, &nat, typelist, atoms->atom[i].type, plist[ftype].interactionTypes, ftype);
346             atoms->atom[i].typeB = search_atomtypes(
347                     this, &nat, typelist, atoms->atom[i].typeB, plist[ftype].interactionTypes, ftype);
348         }
349     }
350
351     for (int i = 0; i < 2; i++)
352     {
353         if (wall_atomtype[i] >= 0)
354         {
355             wall_atomtype[i] = search_atomtypes(
356                     this, &nat, typelist, wall_atomtype[i], plist[ftype].interactionTypes, ftype);
357         }
358     }
359
360     std::vector<AtomTypeData> new_types;
361     /* We now have a list of unique atomtypes in typelist */
362
363     /* Renumber nlist */
364     std::vector<InteractionOfType> nbsnew;
365
366     for (int i = 0; (i < nat); i++)
367     {
368         int mi = typelist[i];
369         for (int j = 0; (j < nat); j++)
370         {
371             int                      mj              = typelist[j];
372             const InteractionOfType& interactionType = plist[ftype].interactionTypes[ntype * mi + mj];
373             nbsnew.emplace_back(interactionType.atoms(),
374                                 interactionType.forceParam(),
375                                 interactionType.interactionTypeName());
376         }
377         new_types.push_back(impl_->types[mi]);
378     }
379
380     mtop->ffparams.atnr = nat;
381
382     impl_->types                  = new_types;
383     plist[ftype].interactionTypes = nbsnew;
384 }
385
386 void PreprocessingAtomTypes::copyTot_atomtypes(t_atomtypes* atomtypes) const
387 {
388     /* Copy the atomtype data to the topology atomtype list */
389     int ntype     = size();
390     atomtypes->nr = ntype;
391     snew(atomtypes->atomnumber, ntype);
392
393     for (int i = 0; i < ntype; i++)
394     {
395         atomtypes->atomnumber[i] = impl_->types[i].atomNumber_;
396     }
397 }