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) 2013,2014,2015,2016,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.
52 #include <unordered_set>
53 #include <sys/types.h>
55 #include "gromacs/fileio/gmxfio.h"
56 #include "gromacs/fileio/warninp.h"
57 #include "gromacs/gmxpreprocess/gmxcpp.h"
58 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
59 #include "gromacs/gmxpreprocess/gpp_bond_atomtype.h"
60 #include "gromacs/gmxpreprocess/gpp_nextnb.h"
61 #include "gromacs/gmxpreprocess/grompp_impl.h"
62 #include "gromacs/gmxpreprocess/readir.h"
63 #include "gromacs/gmxpreprocess/topdirs.h"
64 #include "gromacs/gmxpreprocess/toppush.h"
65 #include "gromacs/gmxpreprocess/topshake.h"
66 #include "gromacs/gmxpreprocess/toputil.h"
67 #include "gromacs/gmxpreprocess/vsite_parm.h"
68 #include "gromacs/math/units.h"
69 #include "gromacs/math/utilities.h"
70 #include "gromacs/mdtypes/inputrec.h"
71 #include "gromacs/mdtypes/md_enums.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/topology/block.h"
74 #include "gromacs/topology/exclusionblocks.h"
75 #include "gromacs/topology/ifunc.h"
76 #include "gromacs/topology/symtab.h"
77 #include "gromacs/topology/topology.h"
78 #include "gromacs/utility/cstringutil.h"
79 #include "gromacs/utility/fatalerror.h"
80 #include "gromacs/utility/futil.h"
81 #include "gromacs/utility/gmxassert.h"
82 #include "gromacs/utility/pleasecite.h"
83 #include "gromacs/utility/smalloc.h"
85 #define OPENDIR '[' /* starting sign for directive */
86 #define CLOSEDIR ']' /* ending sign for directive */
88 static void gen_pairs(const InteractionsOfType& nbs, InteractionsOfType* pairs, real fudge, int comb)
92 int nnn = static_cast<int>(std::sqrt(static_cast<double>(ntp)));
93 GMX_ASSERT(nnn * nnn == ntp,
94 "Number of pairs of generated non-bonded parameters should be a perfect square");
95 int nrfp = NRFP(F_LJ);
96 int nrfpA = interaction_function[F_LJ14].nrfpA;
97 int nrfpB = interaction_function[F_LJ14].nrfpB;
99 if ((nrfp != nrfpA) || (nrfpA != nrfpB))
101 gmx_incons("Number of force parameters in gen_pairs wrong");
104 fprintf(stderr, "Generating 1-4 interactions: fudge = %g\n", fudge);
105 pairs->interactionTypes.clear();
107 std::array<int, 2> atomNumbers;
108 std::array<real, MAXFORCEPARAM> forceParam = { NOTSET };
109 for (const auto& type : nbs.interactionTypes)
111 /* Copy type.atoms */
112 atomNumbers = { i / nnn, i % nnn };
113 /* Copy normal and FEP parameters and multiply by fudge factor */
114 gmx::ArrayRef<const real> existingParam = type.forceParam();
115 GMX_RELEASE_ASSERT(2 * nrfp <= MAXFORCEPARAM,
116 "Can't have more parameters than half of maximum p arameter number");
117 for (int j = 0; j < nrfp; j++)
119 /* If we are using sigma/epsilon values, only the epsilon values
120 * should be scaled, but not sigma.
121 * The sigma values have even indices 0,2, etc.
123 if ((comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) && (j % 2 == 0))
132 forceParam[j] = scaling * existingParam[j];
133 forceParam[nrfp + j] = scaling * existingParam[j];
135 pairs->interactionTypes.emplace_back(InteractionOfType(atomNumbers, forceParam));
140 double check_mol(const gmx_mtop_t* mtop, warninp* wi)
147 /* Check mass and charge */
150 for (const gmx_molblock_t& molb : mtop->molblock)
152 const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
153 for (i = 0; (i < atoms->nr); i++)
155 q += molb.nmol * atoms->atom[i].q;
156 m = atoms->atom[i].m;
157 mB = atoms->atom[i].mB;
158 pt = atoms->atom[i].ptype;
159 /* If the particle is an atom or a nucleus it must have a mass,
160 * else, if it is a shell, a vsite or a bondshell it can have mass zero
162 if (((m <= 0.0) || (mB <= 0.0)) && ((pt == eptAtom) || (pt == eptNucleus)))
164 ri = atoms->atom[i].resind;
165 sprintf(buf, "atom %s (Res %s-%d) has mass %g (state A) / %g (state B)\n",
166 *(atoms->atomname[i]), *(atoms->resinfo[ri].name), atoms->resinfo[ri].nr, m, mB);
167 warning_error(wi, buf);
169 else if (((m != 0) || (mB != 0)) && (pt == eptVSite))
171 ri = atoms->atom[i].resind;
173 "virtual site %s (Res %s-%d) has non-zero mass %g (state A) / %g (state "
175 " Check your topology.\n",
176 *(atoms->atomname[i]), *(atoms->resinfo[ri].name), atoms->resinfo[ri].nr, m, mB);
177 warning_error(wi, buf);
178 /* The following statements make LINCS break! */
179 /* atoms->atom[i].m=0; */
186 /*! \brief Returns the rounded charge of a molecule, when close to integer, otherwise returns the original charge.
188 * The results of this routine are only used for checking and for
189 * printing warning messages. Thus we can assume that charges of molecules
190 * should be integer. If the user wanted non-integer molecular charge,
191 * an undesired warning is printed and the user should use grompp -maxwarn 1.
193 * \param qMol The total, unrounded, charge of the molecule
194 * \param sumAbsQ The sum of absolute values of the charges, used for determining the tolerance for the rounding.
196 static double roundedMoleculeCharge(double qMol, double sumAbsQ)
198 /* We use a tolerance of 1e-6 for inaccuracies beyond the 6th decimal
199 * of the charges for ascii float truncation in the topology files.
200 * Although the summation here uses double precision, the charges
201 * are read and stored in single precision when real=float. This can
202 * lead to rounding errors of half the least significant bit.
203 * Note that, unfortunately, we can not assume addition of random
204 * rounding errors. It is not entirely unlikely that many charges
205 * have a near half-bit rounding error with the same sign.
207 double tolAbs = 1e-6;
208 double tol = std::max(tolAbs, 0.5 * GMX_REAL_EPS * sumAbsQ);
209 double qRound = std::round(qMol);
210 if (std::abs(qMol - qRound) <= tol)
220 static void sum_q(const t_atoms* atoms, int numMols, double* qTotA, double* qTotB)
227 for (int i = 0; i < atoms->nr; i++)
229 qmolA += atoms->atom[i].q;
230 qmolB += atoms->atom[i].qB;
231 sumAbsQA += std::abs(atoms->atom[i].q);
232 sumAbsQB += std::abs(atoms->atom[i].qB);
235 *qTotA += numMols * roundedMoleculeCharge(qmolA, sumAbsQA);
236 *qTotB += numMols * roundedMoleculeCharge(qmolB, sumAbsQB);
239 static void get_nbparm(char* nb_str, char* comb_str, int* nb, int* comb, warninp* wi)
242 char warn_buf[STRLEN];
245 for (i = 1; (i < eNBF_NR); i++)
247 if (gmx_strcasecmp(nb_str, enbf_names[i]) == 0)
254 *nb = strtol(nb_str, nullptr, 10);
256 if ((*nb < 1) || (*nb >= eNBF_NR))
258 sprintf(warn_buf, "Invalid nonbond function selector '%s' using %s", nb_str, enbf_names[1]);
259 warning_error(wi, warn_buf);
263 for (i = 1; (i < eCOMB_NR); i++)
265 if (gmx_strcasecmp(comb_str, ecomb_names[i]) == 0)
272 *comb = strtol(comb_str, nullptr, 10);
274 if ((*comb < 1) || (*comb >= eCOMB_NR))
276 sprintf(warn_buf, "Invalid combination rule selector '%s' using %s", comb_str, ecomb_names[1]);
277 warning_error(wi, warn_buf);
282 static char** cpp_opts(const char* define, const char* include, warninp* wi)
286 const char* cppadds[2];
287 char** cppopts = nullptr;
288 const char* option[2] = { "-D", "-I" };
289 const char* nopt[2] = { "define", "include" };
293 char warn_buf[STRLEN];
296 cppadds[1] = include;
297 for (n = 0; (n < 2); n++)
304 while ((*ptr != '\0') && isspace(*ptr))
309 while ((*rptr != '\0') && !isspace(*rptr))
316 snew(buf, (len + 1));
317 strncpy(buf, ptr, len);
318 if (strstr(ptr, option[n]) != ptr)
320 set_warning_line(wi, "mdp file", -1);
321 sprintf(warn_buf, "Malformed %s option %s", nopt[n], buf);
322 warning(wi, warn_buf);
326 srenew(cppopts, ++ncppopts);
327 cppopts[ncppopts - 1] = gmx_strdup(buf);
335 srenew(cppopts, ++ncppopts);
336 cppopts[ncppopts - 1] = nullptr;
342 static void make_atoms_sys(gmx::ArrayRef<const gmx_molblock_t> molblock,
343 gmx::ArrayRef<const MoleculeInformation> molinfo,
347 atoms->atom = nullptr;
349 for (const gmx_molblock_t& molb : molblock)
351 const t_atoms& mol_atoms = molinfo[molb.type].atoms;
353 srenew(atoms->atom, atoms->nr + molb.nmol * mol_atoms.nr);
355 for (int m = 0; m < molb.nmol; m++)
357 for (int a = 0; a < mol_atoms.nr; a++)
359 atoms->atom[atoms->nr++] = mol_atoms.atom[a];
366 static char** read_topol(const char* infile,
371 PreprocessingAtomTypes* atypes,
372 std::vector<MoleculeInformation>* molinfo,
373 std::unique_ptr<MoleculeInformation>* intermolecular_interactions,
374 gmx::ArrayRef<InteractionsOfType> interactions,
375 int* combination_rule,
379 std::vector<gmx_molblock_t>* molblock,
380 bool* ffParametrizedWithHBondConstraints,
383 bool usingFullRangeElectrostatics,
388 char * pline = nullptr, **title = nullptr;
389 char line[STRLEN], errbuf[256], comb_str[256], nb_str[256];
391 char * dirstr, *dummy2;
392 int nrcopies, nscan, ncombs, ncopy;
393 double fLJ, fQQ, fPOW;
394 MoleculeInformation* mi0 = nullptr;
397 t_nbparam ** nbparam, **pair;
398 real fudgeLJ = -1; /* Multiplication factor to generate 1-4 from LJ */
399 bool bReadDefaults, bReadMolType, bGenPairs, bWarn_copy_A_B;
400 double qt = 0, qBt = 0; /* total charge */
401 int dcatt = -1, nmol_couple;
402 /* File handling variables */
406 char* tmp_line = nullptr;
407 char warn_buf[STRLEN];
408 const char* floating_point_arithmetic_tip =
409 "Total charge should normally be an integer. See\n"
410 "http://www.gromacs.org/Documentation/Floating_Point_Arithmetic\n"
411 "for discussion on how close it should be to an integer.\n";
412 /* We need to open the output file before opening the input file,
413 * because cpp_open_file can change the current working directory.
417 out = gmx_fio_fopen(outfile, "w");
424 /* open input file */
425 auto cpp_opts_return = cpp_opts(define, include, wi);
426 status = cpp_open_file(infile, &handle, cpp_opts_return);
429 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
432 /* some local variables */
433 DS_Init(&DS); /* directive stack */
434 d = Directive::d_invalid; /* first thing should be a directive */
435 nbparam = nullptr; /* The temporary non-bonded matrix */
436 pair = nullptr; /* The temporary pair interaction matrix */
437 std::vector<std::vector<gmx::ExclusionBlock>> exclusionBlocks;
440 *reppow = 12.0; /* Default value for repulsion power */
442 /* Init the number of CMAP torsion angles and grid spacing */
443 interactions[F_CMAP].cmakeGridSpacing = 0;
444 interactions[F_CMAP].cmapAngles = 0;
446 bWarn_copy_A_B = bFEP;
448 PreprocessingBondAtomType bondAtomType;
449 /* parse the actual file */
450 bReadDefaults = FALSE;
452 bReadMolType = FALSE;
457 status = cpp_read_line(&handle, STRLEN, line);
458 done = (status == eCPP_EOF);
461 if (status != eCPP_OK)
463 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
467 fprintf(out, "%s\n", line);
470 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
472 pline = gmx_strdup(line);
474 /* Strip trailing '\' from pline, if it exists */
476 if ((sl > 0) && (pline[sl - 1] == CONTINUE))
481 /* build one long line from several fragments - necessary for CMAP */
482 while (continuing(line))
484 status = cpp_read_line(&handle, STRLEN, line);
485 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
487 /* Since we depend on the '\' being present to continue to read, we copy line
488 * to a tmp string, strip the '\' from that string, and cat it to pline
490 tmp_line = gmx_strdup(line);
492 sl = strlen(tmp_line);
493 if ((sl > 0) && (tmp_line[sl - 1] == CONTINUE))
495 tmp_line[sl - 1] = ' ';
498 done = (status == eCPP_EOF);
501 if (status != eCPP_OK)
503 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
507 fprintf(out, "%s\n", line);
511 srenew(pline, strlen(pline) + strlen(tmp_line) + 1);
512 strcat(pline, tmp_line);
516 /* skip trailing and leading spaces and comment text */
517 strip_comment(pline);
520 /* if there is something left... */
521 if (static_cast<int>(strlen(pline)) > 0)
523 if (pline[0] == OPENDIR)
525 /* A directive on this line: copy the directive
526 * without the brackets into dirstr, then
527 * skip spaces and tabs on either side of directive
529 dirstr = gmx_strdup((pline + 1));
530 if ((dummy2 = strchr(dirstr, CLOSEDIR)) != nullptr)
536 if ((newd = str2dir(dirstr)) == Directive::d_invalid)
538 sprintf(errbuf, "Invalid directive %s", dirstr);
539 warning_error(wi, errbuf);
543 /* Directive found */
544 if (DS_Check_Order(DS, newd))
551 /* we should print here which directives should have
552 been present, and which actually are */
553 gmx_fatal(FARGS, "%s\nInvalid order for directive %s",
554 cpp_error(&handle, eCPP_SYNTAX), dir2str(newd));
555 /* d = Directive::d_invalid; */
558 if (d == Directive::d_intermolecular_interactions)
560 if (*intermolecular_interactions == nullptr)
562 /* We (mis)use the moleculetype processing
563 * to process the intermolecular interactions
564 * by making a "molecule" of the size of the system.
566 *intermolecular_interactions = std::make_unique<MoleculeInformation>();
567 mi0 = intermolecular_interactions->get();
569 make_atoms_sys(*molblock, *molinfo, &mi0->atoms);
575 else if (d != Directive::d_invalid)
577 /* Not a directive, just a plain string
578 * use a gigantic switch to decode,
579 * if there is a valid directive!
583 case Directive::d_defaults:
586 gmx_fatal(FARGS, "%s\nFound a second defaults directive.\n",
587 cpp_error(&handle, eCPP_SYNTAX));
589 bReadDefaults = TRUE;
590 nscan = sscanf(pline, "%s%s%s%lf%lf%lf", nb_str, comb_str, genpairs,
602 get_nbparm(nb_str, comb_str, &nb_funct, combination_rule, wi);
605 bGenPairs = (gmx::equalCaseInsensitive(genpairs, "Y", 1));
606 if (nb_funct != eNBF_LJ && bGenPairs)
609 "Generating pair parameters is only supported "
610 "with LJ non-bonded interactions");
626 nb_funct = ifunc_index(Directive::d_nonbond_params, nb_funct);
629 case Directive::d_atomtypes:
630 push_at(symtab, atypes, &bondAtomType, pline, nb_funct, &nbparam,
631 bGenPairs ? &pair : nullptr, wi);
634 case Directive::d_bondtypes:
635 push_bt(d, interactions, 2, nullptr, &bondAtomType, pline, wi);
637 case Directive::d_constrainttypes:
638 push_bt(d, interactions, 2, nullptr, &bondAtomType, pline, wi);
640 case Directive::d_pairtypes:
643 push_nbt(d, pair, atypes, pline, F_LJ14, wi);
647 push_bt(d, interactions, 2, atypes, nullptr, pline, wi);
650 case Directive::d_angletypes:
651 push_bt(d, interactions, 3, nullptr, &bondAtomType, pline, wi);
653 case Directive::d_dihedraltypes:
654 /* Special routine that can read both 2 and 4 atom dihedral definitions. */
655 push_dihedraltype(d, interactions, &bondAtomType, pline, wi);
658 case Directive::d_nonbond_params:
659 push_nbt(d, nbparam, atypes, pline, nb_funct, wi);
662 case Directive::d_implicit_genborn_params:
663 // Skip this line, so old topologies with
664 // GB parameters can be read.
667 case Directive::d_implicit_surface_params:
668 // Skip this line, so that any topologies
669 // with surface parameters can be read
670 // (even though these were never formally
674 case Directive::d_cmaptypes:
675 push_cmaptype(d, interactions, 5, atypes, &bondAtomType, pline, wi);
678 case Directive::d_moleculetype:
683 if (opts->couple_moltype != nullptr
684 && (opts->couple_lam0 == ecouplamNONE || opts->couple_lam0 == ecouplamQ
685 || opts->couple_lam1 == ecouplamNONE
686 || opts->couple_lam1 == ecouplamQ))
688 dcatt = add_atomtype_decoupled(symtab, atypes, &nbparam,
689 bGenPairs ? &pair : nullptr);
691 ntype = atypes->size();
692 ncombs = (ntype * (ntype + 1)) / 2;
693 generate_nbparams(*combination_rule, nb_funct,
694 &(interactions[nb_funct]), atypes, wi);
695 ncopy = copy_nbparams(nbparam, nb_funct, &(interactions[nb_funct]), ntype);
697 "Generated %d of the %d non-bonded parameter "
699 ncombs - ncopy, ncombs);
700 free_nbparam(nbparam, ntype);
703 gen_pairs((interactions[nb_funct]), &(interactions[F_LJ14]),
704 fudgeLJ, *combination_rule);
705 ncopy = copy_nbparams(pair, nb_funct, &(interactions[F_LJ14]), ntype);
707 "Generated %d of the %d 1-4 parameter combinations\n",
708 ncombs - ncopy, ncombs);
709 free_nbparam(pair, ntype);
711 /* Copy GBSA parameters to atomtype array? */
716 push_molt(symtab, molinfo, pline, wi);
717 exclusionBlocks.emplace_back();
718 mi0 = &molinfo->back();
719 mi0->atoms.haveMass = TRUE;
720 mi0->atoms.haveCharge = TRUE;
721 mi0->atoms.haveType = TRUE;
722 mi0->atoms.haveBState = TRUE;
723 mi0->atoms.havePdbInfo = FALSE;
726 case Directive::d_atoms:
727 push_atom(symtab, &(mi0->atoms), atypes, pline, wi);
730 case Directive::d_pairs:
733 "Need to have a valid MoleculeInformation object to work on");
734 push_bond(d, interactions, mi0->interactions, &(mi0->atoms), atypes,
735 pline, FALSE, bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
737 case Directive::d_pairs_nb:
740 "Need to have a valid MoleculeInformation object to work on");
741 push_bond(d, interactions, mi0->interactions, &(mi0->atoms), atypes,
742 pline, FALSE, FALSE, 1.0, bZero, &bWarn_copy_A_B, wi);
745 case Directive::d_vsites2:
746 case Directive::d_vsites3:
747 case Directive::d_vsites4:
748 case Directive::d_bonds:
749 case Directive::d_angles:
750 case Directive::d_constraints:
751 case Directive::d_settles:
752 case Directive::d_position_restraints:
753 case Directive::d_angle_restraints:
754 case Directive::d_angle_restraints_z:
755 case Directive::d_distance_restraints:
756 case Directive::d_orientation_restraints:
757 case Directive::d_dihedral_restraints:
758 case Directive::d_dihedrals:
759 case Directive::d_polarization:
760 case Directive::d_water_polarization:
761 case Directive::d_thole_polarization:
764 "Need to have a valid MoleculeInformation object to work on");
765 push_bond(d, interactions, mi0->interactions, &(mi0->atoms), atypes,
766 pline, TRUE, bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
768 case Directive::d_cmap:
771 "Need to have a valid MoleculeInformation object to work on");
772 push_cmap(d, interactions, mi0->interactions, &(mi0->atoms), atypes, pline, wi);
775 case Directive::d_vsitesn:
778 "Need to have a valid MoleculeInformation object to work on");
779 push_vsitesn(d, mi0->interactions, &(mi0->atoms), pline, wi);
781 case Directive::d_exclusions:
782 GMX_ASSERT(!exclusionBlocks.empty(),
783 "exclusionBlocks must always be allocated so exclusions can "
785 if (exclusionBlocks.back().empty())
787 GMX_RELEASE_ASSERT(mi0,
788 "Need to have a valid MoleculeInformation "
789 "object to work on");
790 exclusionBlocks.back().resize(mi0->atoms.nr);
792 push_excl(pline, exclusionBlocks.back(), wi);
794 case Directive::d_system:
796 title = put_symtab(symtab, pline);
798 case Directive::d_molecules:
803 push_mol(*molinfo, pline, &whichmol, &nrcopies, wi);
804 mi0 = &((*molinfo)[whichmol]);
805 molblock->resize(molblock->size() + 1);
806 molblock->back().type = whichmol;
807 molblock->back().nmol = nrcopies;
809 bCouple = (opts->couple_moltype != nullptr
810 && (gmx_strcasecmp("system", opts->couple_moltype) == 0
811 || strcmp(*(mi0->name), opts->couple_moltype) == 0));
814 nmol_couple += nrcopies;
817 if (mi0->atoms.nr == 0)
819 gmx_fatal(FARGS, "Molecule type '%s' contains no atoms", *mi0->name);
821 fprintf(stderr, "Excluding %d bonded neighbours molecule type '%s'\n",
822 mi0->nrexcl, *mi0->name);
823 sum_q(&mi0->atoms, nrcopies, &qt, &qBt);
824 if (!mi0->bProcessed)
826 generate_excl(mi0->nrexcl, mi0->atoms.nr, mi0->interactions, &(mi0->excls));
827 gmx::mergeExclusions(&(mi0->excls), exclusionBlocks[whichmol]);
828 make_shake(mi0->interactions, &mi0->atoms, opts->nshake);
832 convert_moltype_couple(mi0, dcatt, *fudgeQQ, opts->couple_lam0,
833 opts->couple_lam1, opts->bCoupleIntra,
834 nb_funct, &(interactions[nb_funct]), wi);
836 stupid_fill_block(&mi0->mols, mi0->atoms.nr, TRUE);
837 mi0->bProcessed = TRUE;
842 fprintf(stderr, "case: %d\n", static_cast<int>(d));
843 gmx_incons("unknown directive");
852 // Check that all strings defined with -D were used when processing topology
853 std::string unusedDefineWarning = checkAndWarnForUnusedDefines(*handle);
854 if (!unusedDefineWarning.empty())
856 warning(wi, unusedDefineWarning);
859 sfree(cpp_opts_return);
866 /* List of GROMACS define names for force fields that have been
867 * parametrized using constraints involving hydrogens only.
869 * We should avoid hardcoded names, but this is hopefully only
870 * needed temparorily for discouraging use of constraints=all-bonds.
872 const std::array<std::string, 3> ffDefines = { "_FF_AMBER", "_FF_CHARMM", "_FF_OPLSAA" };
873 *ffParametrizedWithHBondConstraints = false;
874 for (const std::string& ffDefine : ffDefines)
876 if (cpp_find_define(&handle, ffDefine))
878 *ffParametrizedWithHBondConstraints = true;
882 if (cpp_find_define(&handle, "_FF_GROMOS96") != nullptr)
885 "The GROMOS force fields have been parametrized with a physically incorrect "
886 "multiple-time-stepping scheme for a twin-range cut-off. When used with "
887 "a single-range cut-off (or a correct Trotter multiple-time-stepping scheme), "
888 "physical properties, such as the density, might differ from the intended values. "
889 "Since there are researchers actively working on validating GROMOS with modern "
890 "integrators we have not yet removed the GROMOS force fields, but you should be "
891 "aware of these issues and check if molecules in your system are affected before "
893 "Further information is available at https://redmine.gromacs.org/issues/2884 , "
894 "and a longer explanation of our decision to remove physically incorrect "
896 "can be found at https://doi.org/10.26434/chemrxiv.11474583.v1 .");
901 if (opts->couple_moltype)
903 if (nmol_couple == 0)
905 gmx_fatal(FARGS, "Did not find any molecules of type '%s' for coupling", opts->couple_moltype);
907 fprintf(stderr, "Coupling %d copies of molecule type '%s'\n", nmol_couple, opts->couple_moltype);
910 /* this is not very clean, but fixes core dump on empty system name */
913 title = put_symtab(symtab, "");
918 sprintf(warn_buf, "System has non-zero total charge: %.6f\n%s\n", qt, floating_point_arithmetic_tip);
919 warning_note(wi, warn_buf);
921 if (fabs(qBt) > 1e-4 && !gmx_within_tol(qBt, qt, 1e-6))
923 sprintf(warn_buf, "State B has non-zero total charge: %.6f\n%s\n", qBt,
924 floating_point_arithmetic_tip);
925 warning_note(wi, warn_buf);
927 if (usingFullRangeElectrostatics && (fabs(qt) > 1e-4 || fabs(qBt) > 1e-4))
930 "You are using Ewald electrostatics in a system with net charge. This can lead to "
931 "severe artifacts, such as ions moving into regions with low dielectric, due to "
932 "the uniform background charge. We suggest to neutralize your system with counter "
933 "ions, possibly in combination with a physiological salt concentration.");
934 please_cite(stdout, "Hub2014a");
939 if (*intermolecular_interactions != nullptr)
941 sfree(intermolecular_interactions->get()->atoms.atom);
947 char** do_top(bool bVerbose,
949 const char* topppfile,
953 gmx::ArrayRef<InteractionsOfType> interactions,
954 int* combination_rule,
955 double* repulsion_power,
957 PreprocessingAtomTypes* atypes,
958 std::vector<MoleculeInformation>* molinfo,
959 std::unique_ptr<MoleculeInformation>* intermolecular_interactions,
960 const t_inputrec* ir,
961 std::vector<gmx_molblock_t>* molblock,
962 bool* ffParametrizedWithHBondConstraints,
965 /* Tmpfile might contain a long path */
980 printf("processing topology...\n");
982 title = read_topol(topfile, tmpfile, opts->define, opts->include, symtab, atypes, molinfo,
983 intermolecular_interactions, interactions, combination_rule, repulsion_power,
984 opts, fudgeQQ, molblock, ffParametrizedWithHBondConstraints,
985 ir->efep != efepNO, bZero, EEL_FULL(ir->coulombtype), wi);
987 if ((*combination_rule != eCOMB_GEOMETRIC) && (ir->vdwtype == evdwUSER))
990 "Using sigma/epsilon based combination rules with"
991 " user supplied potential function may produce unwanted"
999 * Generate exclusion lists for QM/MM.
1001 * This routine updates the exclusion lists for QM atoms in order to include all other QM
1002 * atoms of this molecule. Moreover, this routine replaces bonds between QM atoms with
1003 * CONNBOND and, when MiMiC is not used, removes bonded interactions between QM and link atoms.
1004 * Finally, in case if MiMiC QM/MM is used - charges of QM atoms are set to 0
1006 * @param molt molecule type with QM atoms
1007 * @param grpnr group informatio
1008 * @param ir input record
1009 * @param qmmmMode QM/MM mode switch: original/MiMiC
1011 static void generate_qmexcl_moltype(gmx_moltype_t* molt, const unsigned char* grpnr, t_inputrec* ir, GmxQmmmMode qmmmMode)
1013 /* This routine expects molt->ilist to be of size F_NRE and ordered. */
1015 /* generates the exclusions between the individual QM atoms, as
1016 * these interactions should be handled by the QM subroutines and
1017 * not by the gromacs routines
1019 int qm_max = 0, qm_nr = 0, link_nr = 0, link_max = 0;
1020 int * qm_arr = nullptr, *link_arr = nullptr;
1021 bool *bQMMM, *blink;
1023 /* First we search and select the QM atoms in an qm_arr array that
1024 * we use to create the exclusions.
1026 * we take the possibility into account that a user has defined more
1027 * than one QM group:
1029 * for that we also need to do this an ugly work-about just in case
1030 * the QM group contains the entire system...
1033 /* we first search for all the QM atoms and put them in an array
1035 for (int j = 0; j < ir->opts.ngQM; j++)
1037 for (int i = 0; i < molt->atoms.nr; i++)
1039 if (qm_nr >= qm_max)
1042 srenew(qm_arr, qm_max);
1044 if ((grpnr ? grpnr[i] : 0) == j)
1046 qm_arr[qm_nr++] = i;
1047 molt->atoms.atom[i].q = 0.0;
1048 molt->atoms.atom[i].qB = 0.0;
1052 /* bQMMM[..] is an array containin TRUE/FALSE for atoms that are
1053 * QM/not QM. We first set all elements to false. Afterwards we use
1054 * the qm_arr to change the elements corresponding to the QM atoms
1057 snew(bQMMM, molt->atoms.nr);
1058 for (int i = 0; i < molt->atoms.nr; i++)
1062 for (int i = 0; i < qm_nr; i++)
1064 bQMMM[qm_arr[i]] = TRUE;
1067 /* We remove all bonded interactions (i.e. bonds,
1068 * angles, dihedrals, 1-4's), involving the QM atoms. The way they
1069 * are removed is as follows: if the interaction invloves 2 atoms,
1070 * it is removed if both atoms are QMatoms. If it involves 3 atoms,
1071 * it is removed if at least two of the atoms are QM atoms, if the
1072 * interaction involves 4 atoms, it is removed if there are at least
1073 * 2 QM atoms. Since this routine is called once before any forces
1074 * are computed, the top->idef.il[N].iatom[] array (see idef.h) can
1075 * be rewritten at this poitn without any problem. 25-9-2002 */
1077 /* first check whether we already have CONNBONDS.
1078 * Note that if we don't, we don't add a param entry and set ftype=0,
1079 * which is ok, since CONNBONDS does not use parameters.
1081 int ftype_connbond = 0;
1082 int ind_connbond = 0;
1083 if (molt->ilist[F_CONNBONDS].size() != 0)
1085 fprintf(stderr, "nr. of CONNBONDS present already: %d\n", molt->ilist[F_CONNBONDS].size() / 3);
1086 ftype_connbond = molt->ilist[F_CONNBONDS].iatoms[0];
1087 ind_connbond = molt->ilist[F_CONNBONDS].size();
1089 /* now we delete all bonded interactions, except the ones describing
1090 * a chemical bond. These are converted to CONNBONDS
1092 for (int ftype = 0; ftype < F_NRE; ftype++)
1094 if (!(interaction_function[ftype].flags & IF_BOND) || ftype == F_CONNBONDS)
1098 int nratoms = interaction_function[ftype].nratoms;
1100 while (j < molt->ilist[ftype].size())
1106 /* Remove an interaction between two atoms when both are
1107 * in the QM region. Note that we don't have to worry about
1108 * link atoms here, as they won't have 2-atom interactions.
1110 int a1 = molt->ilist[ftype].iatoms[1 + j + 0];
1111 int a2 = molt->ilist[ftype].iatoms[1 + j + 1];
1112 bexcl = (bQMMM[a1] && bQMMM[a2]);
1113 /* A chemical bond between two QM atoms will be copied to
1114 * the F_CONNBONDS list, for reasons mentioned above.
1116 if (bexcl && IS_CHEMBOND(ftype))
1118 InteractionList& ilist = molt->ilist[F_CONNBONDS];
1119 ilist.iatoms.resize(ind_connbond + 3);
1120 ilist.iatoms[ind_connbond++] = ftype_connbond;
1121 ilist.iatoms[ind_connbond++] = a1;
1122 ilist.iatoms[ind_connbond++] = a2;
1127 /* MM interactions have to be excluded if they are included
1128 * in the QM already. Because we use a link atom (H atom)
1129 * when the QM/MM boundary runs through a chemical bond, this
1130 * means that as long as one atom is MM, we still exclude,
1131 * as the interaction is included in the QM via:
1132 * QMatom1-QMatom2-QMatom-3-Linkatom.
1135 for (int jj = j + 1; jj < j + 1 + nratoms; jj++)
1137 if (bQMMM[molt->ilist[ftype].iatoms[jj]])
1143 /* MiMiC treats link atoms as quantum atoms - therefore
1144 * we do not need do additional exclusions here */
1145 if (qmmmMode == GmxQmmmMode::GMX_QMMM_MIMIC)
1147 bexcl = numQmAtoms == nratoms;
1151 bexcl = (numQmAtoms >= nratoms - 1);
1154 if (bexcl && ftype == F_SETTLE)
1157 "Can not apply QM to molecules with SETTLE, replace the moleculetype "
1158 "using QM and SETTLE by one without SETTLE");
1163 /* since the interaction involves QM atoms, these should be
1164 * removed from the MM ilist
1166 InteractionList& ilist = molt->ilist[ftype];
1167 for (int k = j; k < ilist.size() - (nratoms + 1); k++)
1169 ilist.iatoms[k] = ilist.iatoms[k + (nratoms + 1)];
1171 ilist.iatoms.resize(ilist.size() - (nratoms + 1));
1175 j += nratoms + 1; /* the +1 is for the functype */
1179 /* Now, we search for atoms bonded to a QM atom because we also want
1180 * to exclude their nonbonded interactions with the QM atoms. The
1181 * reason for this is that this interaction is accounted for in the
1182 * linkatoms interaction with the QMatoms and would be counted
1185 if (qmmmMode != GmxQmmmMode::GMX_QMMM_MIMIC)
1187 for (int i = 0; i < F_NRE; i++)
1192 while (j < molt->ilist[i].size())
1194 int a1 = molt->ilist[i].iatoms[j + 1];
1195 int a2 = molt->ilist[i].iatoms[j + 2];
1196 if ((bQMMM[a1] && !bQMMM[a2]) || (!bQMMM[a1] && bQMMM[a2]))
1198 if (link_nr >= link_max)
1201 srenew(link_arr, link_max);
1205 link_arr[link_nr++] = a2;
1209 link_arr[link_nr++] = a1;
1217 snew(blink, molt->atoms.nr);
1218 for (int i = 0; i < molt->atoms.nr; i++)
1223 if (qmmmMode != GmxQmmmMode::GMX_QMMM_MIMIC)
1225 for (int i = 0; i < link_nr; i++)
1227 blink[link_arr[i]] = TRUE;
1230 /* creating the exclusion block for the QM atoms. Each QM atom has
1231 * as excluded elements all the other QMatoms (and itself).
1234 qmexcl.nr = molt->atoms.nr;
1235 qmexcl.nra = qm_nr * (qm_nr + link_nr) + link_nr * qm_nr;
1236 snew(qmexcl.index, qmexcl.nr + 1);
1237 snew(qmexcl.a, qmexcl.nra);
1239 for (int i = 0; i < qmexcl.nr; i++)
1241 qmexcl.index[i] = j;
1244 for (int k = 0; k < qm_nr; k++)
1246 qmexcl.a[k + j] = qm_arr[k];
1248 for (int k = 0; k < link_nr; k++)
1250 qmexcl.a[qm_nr + k + j] = link_arr[k];
1252 j += (qm_nr + link_nr);
1256 for (int k = 0; k < qm_nr; k++)
1258 qmexcl.a[k + j] = qm_arr[k];
1263 qmexcl.index[qmexcl.nr] = j;
1265 /* and merging with the exclusions already present in sys.
1268 std::vector<gmx::ExclusionBlock> qmexcl2(molt->atoms.nr);
1269 gmx::blockaToExclusionBlocks(&qmexcl, qmexcl2);
1270 gmx::mergeExclusions(&(molt->excls), qmexcl2);
1272 /* Finally, we also need to get rid of the pair interactions of the
1273 * classical atom bonded to the boundary QM atoms with the QMatoms,
1274 * as this interaction is already accounted for by the QM, so also
1275 * here we run the risk of double counting! We proceed in a similar
1276 * way as we did above for the other bonded interactions: */
1277 for (int i = F_LJ14; i < F_COUL14; i++)
1279 int nratoms = interaction_function[i].nratoms;
1281 while (j < molt->ilist[i].size())
1283 int a1 = molt->ilist[i].iatoms[j + 1];
1284 int a2 = molt->ilist[i].iatoms[j + 2];
1286 ((bQMMM[a1] && bQMMM[a2]) || (blink[a1] && bQMMM[a2]) || (bQMMM[a1] && blink[a2]));
1289 /* since the interaction involves QM atoms, these should be
1290 * removed from the MM ilist
1292 InteractionList& ilist = molt->ilist[i];
1293 for (int k = j; k < ilist.size() - (nratoms + 1); k++)
1295 ilist.iatoms[k] = ilist.iatoms[k + (nratoms + 1)];
1297 ilist.iatoms.resize(ilist.size() - (nratoms + 1));
1301 j += nratoms + 1; /* the +1 is for the functype */
1310 } /* generate_qmexcl */
1312 void generate_qmexcl(gmx_mtop_t* sys, t_inputrec* ir, warninp* wi, GmxQmmmMode qmmmMode)
1314 /* This routine expects molt->molt[m].ilist to be of size F_NRE and ordered.
1317 unsigned char* grpnr;
1318 int mol, nat_mol, nr_mol_with_qm_atoms = 0;
1319 gmx_molblock_t* molb;
1321 int index_offset = 0;
1324 grpnr = sys->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].data();
1326 for (size_t mb = 0; mb < sys->molblock.size(); mb++)
1328 molb = &sys->molblock[mb];
1329 nat_mol = sys->moltype[molb->type].atoms.nr;
1330 for (mol = 0; mol < molb->nmol; mol++)
1333 for (int i = 0; i < nat_mol; i++)
1335 if ((grpnr ? grpnr[i] : 0) < (ir->opts.ngQM))
1344 nr_mol_with_qm_atoms++;
1347 /* We need to split this molblock */
1350 /* Split the molblock at this molecule */
1351 auto pos = sys->molblock.begin() + mb + 1;
1352 sys->molblock.insert(pos, sys->molblock[mb]);
1353 sys->molblock[mb].nmol = mol;
1354 sys->molblock[mb + 1].nmol -= mol;
1356 molb = &sys->molblock[mb];
1360 /* Split the molblock after this molecule */
1361 auto pos = sys->molblock.begin() + mb + 1;
1362 sys->molblock.insert(pos, sys->molblock[mb]);
1363 molb = &sys->molblock[mb];
1364 sys->molblock[mb].nmol = 1;
1365 sys->molblock[mb + 1].nmol -= 1;
1368 /* Create a copy of a moltype for a molecule
1369 * containing QM atoms and append it in the end of the list
1371 std::vector<gmx_moltype_t> temp(sys->moltype.size());
1372 for (size_t i = 0; i < sys->moltype.size(); ++i)
1374 copy_moltype(&sys->moltype[i], &temp[i]);
1376 sys->moltype.resize(sys->moltype.size() + 1);
1377 for (size_t i = 0; i < temp.size(); ++i)
1379 copy_moltype(&temp[i], &sys->moltype[i]);
1381 copy_moltype(&sys->moltype[molb->type], &sys->moltype.back());
1382 /* Copy the exclusions to a new array, since this is the only
1383 * thing that needs to be modified for QMMM.
1385 copy_blocka(&sys->moltype[molb->type].excls, &sys->moltype.back().excls);
1386 /* Set the molecule type for the QMMM molblock */
1387 molb->type = sys->moltype.size() - 1;
1389 generate_qmexcl_moltype(&sys->moltype[molb->type], grpnr, ir, qmmmMode);
1395 index_offset += nat_mol;
1398 if (qmmmMode == GmxQmmmMode::GMX_QMMM_ORIGINAL && nr_mol_with_qm_atoms > 1)
1400 /* generate a warning is there are QM atoms in different
1401 * topologies. In this case it is not possible at this stage to
1402 * mutualy exclude the non-bonded interactions via the
1403 * exclusions (AFAIK). Instead, the user is advised to use the
1404 * energy group exclusions in the mdp file
1407 "\nThe QM subsystem is divided over multiple topologies. "
1408 "The mutual non-bonded interactions cannot be excluded. "
1409 "There are two ways to achieve this:\n\n"
1410 "1) merge the topologies, such that the atoms of the QM "
1411 "subsystem are all present in one single topology file. "
1412 "In this case this warning will dissappear\n\n"
1413 "2) exclude the non-bonded interactions explicitly via the "
1414 "energygrp-excl option in the mdp file. if this is the case "
1415 "this warning may be ignored"