Clean up grompp memory leaks in rotation and pull structures.
[gromacs.git] / src / gromacs / gmxpreprocess / gpp_atomtype.cpp
blob46376534a0941c3609fe59e83bfc014311793cd0
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,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.
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.
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.
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.
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.
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.
38 #include "gmxpre.h"
40 #include "gpp_atomtype.h"
42 #include <climits>
43 #include <cmath>
44 #include <cstring>
46 #include <algorithm>
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"
59 struct AtomTypeData
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)
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_;
82 class PreprocessingAtomTypes::Impl
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;
91 bool PreprocessingAtomTypes::isSet(int nt) const
93 return ((nt >= 0) && (nt < gmx::ssize(*this)));
96 int PreprocessingAtomTypes::atomTypeFromName(const std::string& str) const
98 /* Atom types are always case sensitive */
99 auto found = std::find_if(impl_->types.begin(), impl_->types.end(),
100 [&str](const auto& type) { return str == *type.name_; });
101 if (found == impl_->types.end())
103 return NOTSET;
105 else
107 return std::distance(impl_->types.begin(), found);
111 size_t PreprocessingAtomTypes::size() const
113 return impl_->size();
116 const char* PreprocessingAtomTypes::atomNameFromAtomType(int nt) const
118 return isSet(nt) ? *(impl_->types[nt].name_) : nullptr;
121 real PreprocessingAtomTypes::atomMassFromAtomType(int nt) const
123 return isSet(nt) ? impl_->types[nt].atom_.m : NOTSET;
126 real PreprocessingAtomTypes::atomChargeFromAtomType(int nt) const
128 return isSet(nt) ? impl_->types[nt].atom_.q : NOTSET;
131 int PreprocessingAtomTypes::atomParticleTypeFromAtomType(int nt) const
133 return isSet(nt) ? impl_->types[nt].atom_.ptype : NOTSET;
136 int PreprocessingAtomTypes::bondAtomTypeFromAtomType(int nt) const
138 return isSet(nt) ? impl_->types[nt].bondAtomType_ : NOTSET;
141 int PreprocessingAtomTypes::atomNumberFromAtomType(int nt) const
143 return isSet(nt) ? impl_->types[nt].atomNumber_ : NOTSET;
146 real PreprocessingAtomTypes::atomNonBondedParamFromAtomType(int nt, int param) const
148 if (!isSet(nt))
150 return NOTSET;
152 gmx::ArrayRef<const real> forceParam = impl_->types[nt].nb_.forceParam();
153 if ((param < 0) || (param >= MAXFORCEPARAM))
155 return NOTSET;
157 return forceParam[param];
160 PreprocessingAtomTypes::PreprocessingAtomTypes() : impl_(new Impl) {}
162 PreprocessingAtomTypes::PreprocessingAtomTypes(PreprocessingAtomTypes&& old) noexcept :
163 impl_(std::move(old.impl_))
167 PreprocessingAtomTypes& PreprocessingAtomTypes::operator=(PreprocessingAtomTypes&& old) noexcept
169 impl_ = std::move(old.impl_);
170 return *this;
173 PreprocessingAtomTypes::~PreprocessingAtomTypes() {}
175 int PreprocessingAtomTypes::addType(t_symtab* tab,
176 const t_atom& a,
177 const std::string& name,
178 const InteractionOfType& nb,
179 int bondAtomType,
180 int atomNumber)
182 int position = atomTypeFromName(name);
183 if (position == NOTSET)
185 impl_->types.emplace_back(a, put_symtab(tab, name.c_str()), nb, bondAtomType, atomNumber);
186 return atomTypeFromName(name);
188 else
190 return position;
194 int PreprocessingAtomTypes::setType(int nt,
195 t_symtab* tab,
196 const t_atom& a,
197 const std::string& name,
198 const InteractionOfType& nb,
199 int bondAtomType,
200 int atomNumber)
202 if (!isSet(nt))
204 return NOTSET;
207 impl_->types[nt].atom_ = a;
208 impl_->types[nt].name_ = put_symtab(tab, name.c_str());
209 impl_->types[nt].nb_ = nb;
210 impl_->types[nt].bondAtomType_ = bondAtomType;
211 impl_->types[nt].atomNumber_ = atomNumber;
213 return nt;
216 void PreprocessingAtomTypes::printTypes(FILE* out)
218 fprintf(out, "[ %s ]\n", dir2str(Directive::d_atomtypes));
219 fprintf(out, "; %6s %8s %8s %8s %12s %12s\n", "type", "mass", "charge", "particle", "c6",
220 "c12");
221 for (auto& entry : impl_->types)
223 fprintf(out, "%8s %8.3f %8.3f %8s %12e %12e\n", *(entry.name_), entry.atom_.m,
224 entry.atom_.q, "A", entry.nb_.c0(), entry.nb_.c1());
227 fprintf(out, "\n");
230 static int search_atomtypes(const PreprocessingAtomTypes* ga,
231 int* n,
232 gmx::ArrayRef<int> typelist,
233 int thistype,
234 gmx::ArrayRef<const InteractionOfType> interactionTypes,
235 int ftype)
237 int nn = *n;
238 int nrfp = NRFP(ftype);
239 int ntype = ga->size();
241 int i;
242 for (i = 0; (i < nn); i++)
244 if (typelist[i] == thistype)
246 /* This type number has already been added */
247 break;
250 /* Otherwise, check if the parameters are identical to any previously added type */
252 bool bFound = true;
253 for (int j = 0; j < ntype && bFound; j++)
255 /* Check nonbonded parameters */
256 gmx::ArrayRef<const real> forceParam1 =
257 interactionTypes[ntype * typelist[i] + j].forceParam();
258 gmx::ArrayRef<const real> forceParam2 = interactionTypes[ntype * thistype + j].forceParam();
259 for (int k = 0; (k < nrfp) && bFound; k++)
261 bFound = forceParam1[k] == forceParam2[k];
264 /* Check atomnumber */
265 int tli = typelist[i];
266 bFound = bFound && (ga->atomNumberFromAtomType(tli) == ga->atomNumberFromAtomType(thistype));
268 if (bFound)
270 break;
274 if (i == nn)
276 if (nn == ntype)
278 gmx_fatal(FARGS, "Atomtype horror n = %d, %s, %d", nn, __FILE__, __LINE__);
280 typelist[nn] = thistype;
281 nn++;
283 *n = nn;
285 return i;
288 void PreprocessingAtomTypes::renumberTypes(gmx::ArrayRef<InteractionsOfType> plist,
289 gmx_mtop_t* mtop,
290 int* wall_atomtype,
291 bool bVerbose)
293 int nat, ftype, ntype;
295 ntype = size();
296 std::vector<int> typelist(ntype);
298 if (bVerbose)
300 fprintf(stderr, "renumbering atomtypes...\n");
303 /* Since the bonded interactions have been assigned now,
304 * we want to reduce the number of atom types by merging
305 * ones with identical nonbonded interactions, in addition
306 * to removing unused ones.
308 * With QM/MM we also check that the atom numbers match
311 /* Get nonbonded interaction type */
312 if (plist[F_LJ].size() > 0)
314 ftype = F_LJ;
316 else
318 ftype = F_BHAM;
321 /* Renumber atomtypes by first making a list of which ones are actually used.
322 * We provide the list of nonbonded parameters so search_atomtypes
323 * can determine if two types should be merged.
325 nat = 0;
326 for (const gmx_moltype_t& moltype : mtop->moltype)
328 const t_atoms* atoms = &moltype.atoms;
329 for (int i = 0; (i < atoms->nr); i++)
331 atoms->atom[i].type = search_atomtypes(this, &nat, typelist, atoms->atom[i].type,
332 plist[ftype].interactionTypes, ftype);
333 atoms->atom[i].typeB = search_atomtypes(this, &nat, typelist, atoms->atom[i].typeB,
334 plist[ftype].interactionTypes, ftype);
338 for (int i = 0; i < 2; i++)
340 if (wall_atomtype[i] >= 0)
342 wall_atomtype[i] = search_atomtypes(this, &nat, typelist, wall_atomtype[i],
343 plist[ftype].interactionTypes, ftype);
347 std::vector<AtomTypeData> new_types;
348 /* We now have a list of unique atomtypes in typelist */
350 /* Renumber nlist */
351 std::vector<InteractionOfType> nbsnew;
353 for (int i = 0; (i < nat); i++)
355 int mi = typelist[i];
356 for (int j = 0; (j < nat); j++)
358 int mj = typelist[j];
359 const InteractionOfType& interactionType = plist[ftype].interactionTypes[ntype * mi + mj];
360 nbsnew.emplace_back(interactionType.atoms(), interactionType.forceParam(),
361 interactionType.interactionTypeName());
363 new_types.push_back(impl_->types[mi]);
366 mtop->ffparams.atnr = nat;
368 impl_->types = new_types;
369 plist[ftype].interactionTypes = nbsnew;
372 void PreprocessingAtomTypes::copyTot_atomtypes(t_atomtypes* atomtypes) const
374 /* Copy the atomtype data to the topology atomtype list */
375 int ntype = size();
376 atomtypes->nr = ntype;
377 snew(atomtypes->atomnumber, ntype);
379 for (int i = 0; i < ntype; i++)
381 atomtypes->atomnumber[i] = impl_->types[i].atomNumber_;