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,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.
47 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
48 #include "gromacs/gmxpreprocess/grompp_impl.h"
49 #include "gromacs/gmxpreprocess/notset.h"
50 #include "gromacs/gmxpreprocess/topdirs.h"
51 #include "gromacs/topology/block.h"
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/topology/symtab.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/gmxassert.h"
56 #include "gromacs/utility/smalloc.h"
60 void add_param_to_list(InteractionsOfType
*list
, const InteractionOfType
&b
)
62 list
->interactionTypes
.emplace_back(b
);
65 /* PRINTING STRUCTURES */
67 static void print_bt(FILE *out
, Directive d
, PreprocessingAtomTypes
*at
,
68 int ftype
, int fsubtype
, gmx::ArrayRef
<const InteractionsOfType
> plist
,
71 /* This dihp is a DIRTY patch because the dih-types do not use
72 * all four atoms to determine the type.
74 const int dihp
[2][2] = { { 1, 2 }, { 0, 3 } };
76 bool bDih
= false, bSwapParity
;
78 const InteractionsOfType
*bt
= &(plist
[ftype
]);
106 case F_CROSS_BOND_ANGLES
:
109 case F_CROSS_BOND_BONDS
:
159 fprintf(out
, "[ %s ]\n", dir2str(d
));
163 fprintf (out
, "%3s %4s", "ai", "aj");
164 for (int j
= 2; (j
< nral
); j
++)
166 fprintf (out
, " %3c%c", 'a', 'i'+j
);
171 for (int j
= 0; (j
< 2); j
++)
173 fprintf (out
, "%3c%c", 'a', 'i'+dihp
[f
][j
]);
177 fprintf (out
, " funct");
178 for (int j
= 0; (j
< nrfp
); j
++)
180 fprintf (out
, " %12c%1d", 'c', j
);
184 /* print bondtypes */
185 for (const auto &parm
: bt
->interactionTypes
)
187 bSwapParity
= (parm
.c0() == NOTSET
) && (parm
.c1() == -1);
188 gmx::ArrayRef
<const int> atoms
= parm
.atoms();
191 for (int j
= 0; (j
< nral
); j
++)
193 fprintf (out
, "%5s ", at
->atomNameFromAtomType(atoms
[j
]));
198 for (int j
= 0; (j
< 2); j
++)
200 fprintf (out
, "%5s ", at
->atomNameFromAtomType(atoms
[dihp
[f
][j
]]));
203 fprintf (out
, "%5d ", bSwapParity
? -f
-1 : f
+1);
205 if (!parm
.interactionTypeName().empty())
207 fprintf(out
, " %s", parm
.interactionTypeName().c_str());
211 gmx::ArrayRef
<const real
> forceParam
= parm
.forceParam();
212 for (int j
= 0; (j
< nrfp
) && (forceParam
[j
] != NOTSET
); j
++)
214 fprintf (out
, "%13.6e ", forceParam
[j
]);
224 void print_blocka(FILE *out
, const char *szName
,
225 const char *szIndex
, const char *szA
,
230 fprintf (out
, "; %s\n", szName
);
231 fprintf (out
, "; %4s %s\n", szIndex
, szA
);
232 for (i
= 0; (i
< block
->nr
); i
++)
234 for (i
= 0; (i
< block
->nr
); i
++)
236 fprintf (out
, "%6d", i
+1);
237 for (j
= block
->index
[i
]; (j
< (block
->index
[i
+1])); j
++)
239 fprintf (out
, "%5d", block
->a
[j
]+1);
247 void print_excl(FILE *out
, int natoms
, t_excls excls
[])
254 for (i
= 0; i
< natoms
&& !have_excl
; i
++)
256 have_excl
= (excls
[i
].nr
> 0);
261 fprintf (out
, "[ %s ]\n", dir2str(Directive::d_exclusions
));
262 fprintf (out
, "; %4s %s\n", "i", "excluded from i");
263 for (i
= 0; i
< natoms
; i
++)
267 fprintf (out
, "%6d ", i
+1);
268 for (j
= 0; j
< excls
[i
].nr
; j
++)
270 fprintf (out
, " %5d", excls
[i
].e
[j
]+1);
280 static double get_residue_charge(const t_atoms
*atoms
, int at
)
285 ri
= atoms
->atom
[at
].resind
;
287 while (at
< atoms
->nr
&& atoms
->atom
[at
].resind
== ri
)
289 q
+= atoms
->atom
[at
].q
;
296 void print_atoms(FILE *out
, PreprocessingAtomTypes
*atype
, t_atoms
*at
, int *cgnr
,
302 const char *tpnmA
, *tpnmB
;
305 as
= dir2str(Directive::d_atoms
);
306 fprintf(out
, "[ %s ]\n", as
);
307 fprintf(out
, "; %4s %10s %6s %7s%6s %6s %10s %10s %6s %10s %10s\n",
308 "nr", "type", "resnr", "residue", "atom", "cgnr", "charge", "mass", "typeB", "chargeB", "massB");
314 /* if the information is present... */
315 for (i
= 0; (i
< at
->nr
); i
++)
317 ri
= at
->atom
[i
].resind
;
318 if ((i
== 0 || ri
!= at
->atom
[i
-1].resind
) &&
319 at
->resinfo
[ri
].rtp
!= nullptr)
321 qres
= get_residue_charge(at
, i
);
322 fprintf(out
, "; residue %3d %-3s rtp %-4s q ",
324 *at
->resinfo
[ri
].name
,
325 *at
->resinfo
[ri
].rtp
);
326 if (fabs(qres
) < 0.001)
328 fprintf(out
, " %s", "0.0");
332 fprintf(out
, "%+3.1f", qres
);
336 tpA
= at
->atom
[i
].type
;
337 if ((tpnmA
= atype
->atomNameFromAtomType(tpA
)) == nullptr)
339 gmx_fatal(FARGS
, "tpA = %d, i= %d in print_atoms", tpA
, i
);
342 /* This is true by construction, but static analysers don't know */
343 GMX_ASSERT(!bRTPresname
|| at
->resinfo
[at
->atom
[i
].resind
].rtp
, "-rtpres did not have residue name available");
344 fprintf(out
, "%6d %10s %6d%c %5s %6s %6d %10g %10g",
349 *(at
->resinfo
[at
->atom
[i
].resind
].rtp
) :
350 *(at
->resinfo
[at
->atom
[i
].resind
].name
),
351 *(at
->atomname
[i
]), cgnr
[i
],
352 at
->atom
[i
].q
, at
->atom
[i
].m
);
353 if (PERTURBED(at
->atom
[i
]))
355 tpB
= at
->atom
[i
].typeB
;
356 if ((tpnmB
= atype
->atomNameFromAtomType(tpB
)) == nullptr)
358 gmx_fatal(FARGS
, "tpB = %d, i= %d in print_atoms", tpB
, i
);
360 fprintf(out
, " %6s %10g %10g",
361 tpnmB
, at
->atom
[i
].qB
, at
->atom
[i
].mB
);
363 // Accumulate the total charge to help troubleshoot issues.
364 qtot
+= static_cast<double>(at
->atom
[i
].q
);
365 // Round it to zero if it is close to zero, because
366 // printing -9.34e-5 confuses users.
367 if (fabs(qtot
) < 0.0001)
371 // Write the total charge for the last atom of the system
372 // and/or residue, because generally that's where it is
373 // expected to be an integer.
374 if (i
== at
->nr
-1 || ri
!= at
->atom
[i
+1].resind
)
376 fprintf(out
, " ; qtot %.4g\n", qtot
);
388 void print_bondeds(FILE *out
, int natoms
, Directive d
,
389 int ftype
, int fsubtype
, gmx::ArrayRef
<const InteractionsOfType
> plist
)
394 PreprocessingAtomTypes atype
;
397 for (int i
= 0; (i
< natoms
); i
++)
400 sprintf(buf
, "%4d", (i
+1));
401 atype
.addType(&stab
, *a
, buf
, InteractionOfType({}, {}), 0, 0);
403 print_bt(out
, d
, &atype
, ftype
, fsubtype
, plist
, TRUE
);