Clean up grompp memory leaks in rotation and pull structures.
[gromacs.git] / src / gromacs / gmxpreprocess / toputil.cpp
blob10ebbfce3e417c6d23309bca8db42c32dc63d2c6
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) 2012,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 "toputil.h"
42 #include <climits>
43 #include <cmath>
44 #include <cstring>
46 #include <algorithm>
48 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
49 #include "gromacs/gmxpreprocess/grompp_impl.h"
50 #include "gromacs/gmxpreprocess/notset.h"
51 #include "gromacs/gmxpreprocess/topdirs.h"
52 #include "gromacs/topology/block.h"
53 #include "gromacs/topology/ifunc.h"
54 #include "gromacs/topology/symtab.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/gmxassert.h"
57 #include "gromacs/utility/smalloc.h"
59 /* UTILITIES */
61 void add_param_to_list(InteractionsOfType* list, const InteractionOfType& b)
63 list->interactionTypes.emplace_back(b);
66 /* PRINTING STRUCTURES */
68 static void print_bt(FILE* out,
69 Directive d,
70 PreprocessingAtomTypes* at,
71 int ftype,
72 int fsubtype,
73 gmx::ArrayRef<const InteractionsOfType> plist,
74 bool bFullDih)
76 /* This dihp is a DIRTY patch because the dih-types do not use
77 * all four atoms to determine the type.
79 const int dihp[2][2] = { { 1, 2 }, { 0, 3 } };
80 int nral, nrfp;
81 bool bDih = false, bSwapParity;
83 const InteractionsOfType* bt = &(plist[ftype]);
85 if (bt->size() == 0)
87 return;
90 int f = 0;
91 switch (ftype)
93 case F_G96ANGLES: // Intended to fall through
94 case F_G96BONDS: f = 1; break;
95 case F_MORSE: f = 2; break;
96 case F_CUBICBONDS: f = 3; break;
97 case F_CONNBONDS: f = 4; break;
98 case F_HARMONIC: f = 5; break;
99 case F_CROSS_BOND_ANGLES: f = 2; break;
100 case F_CROSS_BOND_BONDS: f = 3; break;
101 case F_UREY_BRADLEY: f = 4; break;
102 case F_PDIHS: // Intended to fall through
103 case F_RBDIHS: // Intended to fall through
104 case F_FOURDIHS: bDih = TRUE; break;
105 case F_IDIHS:
106 f = 1;
107 bDih = TRUE;
108 break;
109 case F_CONSTRNC: // Intended to fall through
110 case F_VSITE3FD: f = 1; break;
111 case F_VSITE3FAD: f = 2; break;
112 case F_VSITE3OUT: f = 3; break;
113 case F_VSITE4FDN: // Intended to fall through
114 case F_CMAP: f = 1; break;
116 default: bDih = FALSE;
118 if (bFullDih)
120 bDih = FALSE;
122 if (fsubtype)
124 f = fsubtype - 1;
127 nral = NRAL(ftype);
128 nrfp = NRFP(ftype);
130 /* header */
131 fprintf(out, "[ %s ]\n", dir2str(d));
132 fprintf(out, "; ");
133 if (!bDih)
135 fprintf(out, "%3s %4s", "ai", "aj");
136 for (int j = 2; (j < nral); j++)
138 fprintf(out, " %3c%c", 'a', 'i' + j);
141 else
143 for (int j = 0; (j < 2); j++)
145 fprintf(out, "%3c%c", 'a', 'i' + dihp[f][j]);
149 fprintf(out, " funct");
150 for (int j = 0; (j < nrfp); j++)
152 fprintf(out, " %12c%1d", 'c', j);
154 fprintf(out, "\n");
156 /* print bondtypes */
157 for (const auto& parm : bt->interactionTypes)
159 bSwapParity = (parm.c0() == NOTSET) && (parm.c1() == -1);
160 gmx::ArrayRef<const int> atoms = parm.atoms();
161 if (!bDih)
163 for (int j = 0; (j < nral); j++)
165 fprintf(out, "%5s ", at->atomNameFromAtomType(atoms[j]));
168 else
170 for (int j = 0; (j < 2); j++)
172 fprintf(out, "%5s ", at->atomNameFromAtomType(atoms[dihp[f][j]]));
175 fprintf(out, "%5d ", bSwapParity ? -f - 1 : f + 1);
177 if (!parm.interactionTypeName().empty())
179 fprintf(out, " %s", parm.interactionTypeName().c_str());
181 else
183 gmx::ArrayRef<const real> forceParam = parm.forceParam();
184 for (int j = 0; (j < nrfp) && (forceParam[j] != NOTSET); j++)
186 fprintf(out, "%13.6e ", forceParam[j]);
190 fprintf(out, "\n");
192 fprintf(out, "\n");
193 fflush(out);
196 void print_excl(FILE* out, int natoms, t_excls excls[])
198 int i;
199 bool have_excl;
200 int j;
202 have_excl = FALSE;
203 for (i = 0; i < natoms && !have_excl; i++)
205 have_excl = (excls[i].nr > 0);
208 if (have_excl)
210 fprintf(out, "[ %s ]\n", dir2str(Directive::d_exclusions));
211 fprintf(out, "; %4s %s\n", "i", "excluded from i");
212 for (i = 0; i < natoms; i++)
214 if (excls[i].nr > 0)
216 fprintf(out, "%6d ", i + 1);
217 for (j = 0; j < excls[i].nr; j++)
219 fprintf(out, " %5d", excls[i].e[j] + 1);
221 fprintf(out, "\n");
224 fprintf(out, "\n");
225 fflush(out);
229 static double get_residue_charge(const t_atoms* atoms, int at)
231 int ri;
232 double q;
234 ri = atoms->atom[at].resind;
235 q = 0;
236 while (at < atoms->nr && atoms->atom[at].resind == ri)
238 q += atoms->atom[at].q;
239 at++;
242 return q;
245 void print_atoms(FILE* out, PreprocessingAtomTypes* atype, t_atoms* at, int* cgnr, bool bRTPresname)
247 int i, ri;
248 int tpA, tpB;
249 const char* as;
250 const char *tpnmA, *tpnmB;
251 double qres, qtot;
253 as = dir2str(Directive::d_atoms);
254 fprintf(out, "[ %s ]\n", as);
255 fprintf(out, "; %4s %10s %6s %7s%6s %6s %10s %10s %6s %10s %10s\n", "nr", "type", "resnr",
256 "residue", "atom", "cgnr", "charge", "mass", "typeB", "chargeB", "massB");
258 qtot = 0;
260 if (at->nres)
262 /* if the information is present... */
263 for (i = 0; (i < at->nr); i++)
265 ri = at->atom[i].resind;
266 if ((i == 0 || ri != at->atom[i - 1].resind) && at->resinfo[ri].rtp != nullptr)
268 qres = get_residue_charge(at, i);
269 fprintf(out, "; residue %3d %-3s rtp %-4s q ", at->resinfo[ri].nr,
270 *at->resinfo[ri].name, *at->resinfo[ri].rtp);
271 if (fabs(qres) < 0.001)
273 fprintf(out, " %s", "0.0");
275 else
277 fprintf(out, "%+3.1f", qres);
279 fprintf(out, "\n");
281 tpA = at->atom[i].type;
282 if ((tpnmA = atype->atomNameFromAtomType(tpA)) == nullptr)
284 gmx_fatal(FARGS, "tpA = %d, i= %d in print_atoms", tpA, i);
287 /* This is true by construction, but static analysers don't know */
288 GMX_ASSERT(!bRTPresname || at->resinfo[at->atom[i].resind].rtp,
289 "-rtpres did not have residue name available");
290 fprintf(out, "%6d %10s %6d%c %5s %6s %6d %10g %10g", i + 1, tpnmA, at->resinfo[ri].nr,
291 at->resinfo[ri].ic,
292 bRTPresname ? *(at->resinfo[at->atom[i].resind].rtp)
293 : *(at->resinfo[at->atom[i].resind].name),
294 *(at->atomname[i]), cgnr[i], at->atom[i].q, at->atom[i].m);
295 if (PERTURBED(at->atom[i]))
297 tpB = at->atom[i].typeB;
298 if ((tpnmB = atype->atomNameFromAtomType(tpB)) == nullptr)
300 gmx_fatal(FARGS, "tpB = %d, i= %d in print_atoms", tpB, i);
302 fprintf(out, " %6s %10g %10g", tpnmB, at->atom[i].qB, at->atom[i].mB);
304 // Accumulate the total charge to help troubleshoot issues.
305 qtot += static_cast<double>(at->atom[i].q);
306 // Round it to zero if it is close to zero, because
307 // printing -9.34e-5 confuses users.
308 if (fabs(qtot) < 0.0001)
310 qtot = 0;
312 // Write the total charge for the last atom of the system
313 // and/or residue, because generally that's where it is
314 // expected to be an integer.
315 if (i == at->nr - 1 || ri != at->atom[i + 1].resind)
317 fprintf(out, " ; qtot %.4g\n", qtot);
319 else
321 fputs("\n", out);
325 fprintf(out, "\n");
326 fflush(out);
329 void print_bondeds(FILE* out, int natoms, Directive d, int ftype, int fsubtype, gmx::ArrayRef<const InteractionsOfType> plist)
331 t_symtab stab;
332 t_atom* a;
334 PreprocessingAtomTypes atype;
335 snew(a, 1);
336 open_symtab(&stab);
337 for (int i = 0; (i < natoms); i++)
339 char buf[12];
340 sprintf(buf, "%4d", (i + 1));
341 atype.addType(&stab, *a, buf, InteractionOfType({}, {}), 0, 0);
343 print_bt(out, d, &atype, ftype, fsubtype, plist, TRUE);
345 done_symtab(&stab);
346 sfree(a);