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 by the GROMACS development team.
7 * Copyright (c) 2018,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.
53 #include <unordered_set>
54 #include <sys/types.h>
56 #include "gromacs/fileio/gmxfio.h"
57 #include "gromacs/fileio/warninp.h"
58 #include "gromacs/gmxpreprocess/gmxcpp.h"
59 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
60 #include "gromacs/gmxpreprocess/gpp_bond_atomtype.h"
61 #include "gromacs/gmxpreprocess/gpp_nextnb.h"
62 #include "gromacs/gmxpreprocess/grompp_impl.h"
63 #include "gromacs/gmxpreprocess/readir.h"
64 #include "gromacs/gmxpreprocess/topdirs.h"
65 #include "gromacs/gmxpreprocess/toppush.h"
66 #include "gromacs/gmxpreprocess/topshake.h"
67 #include "gromacs/gmxpreprocess/toputil.h"
68 #include "gromacs/gmxpreprocess/vsite_parm.h"
69 #include "gromacs/math/units.h"
70 #include "gromacs/math/utilities.h"
71 #include "gromacs/mdtypes/inputrec.h"
72 #include "gromacs/mdtypes/md_enums.h"
73 #include "gromacs/pbcutil/pbc.h"
74 #include "gromacs/topology/block.h"
75 #include "gromacs/topology/exclusionblocks.h"
76 #include "gromacs/topology/ifunc.h"
77 #include "gromacs/topology/symtab.h"
78 #include "gromacs/topology/topology.h"
79 #include "gromacs/utility/cstringutil.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/futil.h"
82 #include "gromacs/utility/gmxassert.h"
83 #include "gromacs/utility/logger.h"
84 #include "gromacs/utility/pleasecite.h"
85 #include "gromacs/utility/smalloc.h"
87 #define OPENDIR '[' /* starting sign for directive */
88 #define CLOSEDIR ']' /* ending sign for directive */
90 static void gen_pairs(const InteractionsOfType& nbs, InteractionsOfType* pairs, real fudge, int comb)
94 int nnn = static_cast<int>(std::sqrt(static_cast<double>(ntp)));
95 GMX_ASSERT(nnn * nnn == ntp,
96 "Number of pairs of generated non-bonded parameters should be a perfect square");
97 int nrfp = NRFP(F_LJ);
98 int nrfpA = interaction_function[F_LJ14].nrfpA;
99 int nrfpB = interaction_function[F_LJ14].nrfpB;
101 if ((nrfp != nrfpA) || (nrfpA != nrfpB))
103 gmx_incons("Number of force parameters in gen_pairs wrong");
106 fprintf(stderr, "Generating 1-4 interactions: fudge = %g\n", fudge);
107 pairs->interactionTypes.clear();
109 std::array<int, 2> atomNumbers;
110 std::array<real, MAXFORCEPARAM> forceParam = { NOTSET };
111 for (const auto& type : nbs.interactionTypes)
113 /* Copy type.atoms */
114 atomNumbers = { i / nnn, i % nnn };
115 /* Copy normal and FEP parameters and multiply by fudge factor */
116 gmx::ArrayRef<const real> existingParam = type.forceParam();
117 GMX_RELEASE_ASSERT(2 * nrfp <= MAXFORCEPARAM,
118 "Can't have more parameters than half of maximum p arameter number");
119 for (int j = 0; j < nrfp; j++)
121 /* If we are using sigma/epsilon values, only the epsilon values
122 * should be scaled, but not sigma.
123 * The sigma values have even indices 0,2, etc.
125 if ((comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) && (j % 2 == 0))
134 forceParam[j] = scaling * existingParam[j];
135 forceParam[nrfp + j] = scaling * existingParam[j];
137 pairs->interactionTypes.emplace_back(InteractionOfType(atomNumbers, forceParam));
142 double check_mol(const gmx_mtop_t* mtop, warninp* wi)
149 /* Check mass and charge */
152 for (const gmx_molblock_t& molb : mtop->molblock)
154 const t_atoms* atoms = &mtop->moltype[molb.type].atoms;
155 for (i = 0; (i < atoms->nr); i++)
157 q += molb.nmol * atoms->atom[i].q;
158 m = atoms->atom[i].m;
159 mB = atoms->atom[i].mB;
160 pt = atoms->atom[i].ptype;
161 /* If the particle is an atom or a nucleus it must have a mass,
162 * else, if it is a shell, a vsite or a bondshell it can have mass zero
164 if (((m <= 0.0) || (mB <= 0.0)) && ((pt == eptAtom) || (pt == eptNucleus)))
166 ri = atoms->atom[i].resind;
167 sprintf(buf, "atom %s (Res %s-%d) has mass %g (state A) / %g (state B)\n",
168 *(atoms->atomname[i]), *(atoms->resinfo[ri].name), atoms->resinfo[ri].nr, m, mB);
169 warning_error(wi, buf);
171 else if (((m != 0) || (mB != 0)) && (pt == eptVSite))
173 ri = atoms->atom[i].resind;
175 "virtual site %s (Res %s-%d) has non-zero mass %g (state A) / %g (state "
177 " Check your topology.\n",
178 *(atoms->atomname[i]), *(atoms->resinfo[ri].name), atoms->resinfo[ri].nr, m, mB);
179 warning_error(wi, buf);
180 /* The following statements make LINCS break! */
181 /* atoms->atom[i].m=0; */
188 /*! \brief Returns the rounded charge of a molecule, when close to integer, otherwise returns the original charge.
190 * The results of this routine are only used for checking and for
191 * printing warning messages. Thus we can assume that charges of molecules
192 * should be integer. If the user wanted non-integer molecular charge,
193 * an undesired warning is printed and the user should use grompp -maxwarn 1.
195 * \param qMol The total, unrounded, charge of the molecule
196 * \param sumAbsQ The sum of absolute values of the charges, used for determining the tolerance for the rounding.
198 static double roundedMoleculeCharge(double qMol, double sumAbsQ)
200 /* We use a tolerance of 1e-6 for inaccuracies beyond the 6th decimal
201 * of the charges for ascii float truncation in the topology files.
202 * Although the summation here uses double precision, the charges
203 * are read and stored in single precision when real=float. This can
204 * lead to rounding errors of half the least significant bit.
205 * Note that, unfortunately, we can not assume addition of random
206 * rounding errors. It is not entirely unlikely that many charges
207 * have a near half-bit rounding error with the same sign.
209 double tolAbs = 1e-6;
210 double tol = std::max(tolAbs, 0.5 * GMX_REAL_EPS * sumAbsQ);
211 double qRound = std::round(qMol);
212 if (std::abs(qMol - qRound) <= tol)
222 static void sum_q(const t_atoms* atoms, int numMols, double* qTotA, double* qTotB)
229 for (int i = 0; i < atoms->nr; i++)
231 qmolA += atoms->atom[i].q;
232 qmolB += atoms->atom[i].qB;
233 sumAbsQA += std::abs(atoms->atom[i].q);
234 sumAbsQB += std::abs(atoms->atom[i].qB);
237 *qTotA += numMols * roundedMoleculeCharge(qmolA, sumAbsQA);
238 *qTotB += numMols * roundedMoleculeCharge(qmolB, sumAbsQB);
241 static void get_nbparm(char* nb_str, char* comb_str, int* nb, int* comb, warninp* wi)
244 char warn_buf[STRLEN];
247 for (i = 1; (i < eNBF_NR); i++)
249 if (gmx_strcasecmp(nb_str, enbf_names[i]) == 0)
256 *nb = strtol(nb_str, nullptr, 10);
258 if ((*nb < 1) || (*nb >= eNBF_NR))
260 sprintf(warn_buf, "Invalid nonbond function selector '%s' using %s", nb_str, enbf_names[1]);
261 warning_error(wi, warn_buf);
265 for (i = 1; (i < eCOMB_NR); i++)
267 if (gmx_strcasecmp(comb_str, ecomb_names[i]) == 0)
274 *comb = strtol(comb_str, nullptr, 10);
276 if ((*comb < 1) || (*comb >= eCOMB_NR))
278 sprintf(warn_buf, "Invalid combination rule selector '%s' using %s", comb_str, ecomb_names[1]);
279 warning_error(wi, warn_buf);
284 static char** cpp_opts(const char* define, const char* include, warninp* wi)
288 const char* cppadds[2];
289 char** cppopts = nullptr;
290 const char* option[2] = { "-D", "-I" };
291 const char* nopt[2] = { "define", "include" };
295 char warn_buf[STRLEN];
298 cppadds[1] = include;
299 for (n = 0; (n < 2); n++)
306 while ((*ptr != '\0') && isspace(*ptr))
311 while ((*rptr != '\0') && !isspace(*rptr))
318 snew(buf, (len + 1));
319 strncpy(buf, ptr, len);
320 if (strstr(ptr, option[n]) != ptr)
322 set_warning_line(wi, "mdp file", -1);
323 sprintf(warn_buf, "Malformed %s option %s", nopt[n], buf);
324 warning(wi, warn_buf);
328 srenew(cppopts, ++ncppopts);
329 cppopts[ncppopts - 1] = gmx_strdup(buf);
337 srenew(cppopts, ++ncppopts);
338 cppopts[ncppopts - 1] = nullptr;
344 static void make_atoms_sys(gmx::ArrayRef<const gmx_molblock_t> molblock,
345 gmx::ArrayRef<const MoleculeInformation> molinfo,
349 atoms->atom = nullptr;
351 for (const gmx_molblock_t& molb : molblock)
353 const t_atoms& mol_atoms = molinfo[molb.type].atoms;
355 srenew(atoms->atom, atoms->nr + molb.nmol * mol_atoms.nr);
357 for (int m = 0; m < molb.nmol; m++)
359 for (int a = 0; a < mol_atoms.nr; a++)
361 atoms->atom[atoms->nr++] = mol_atoms.atom[a];
368 static char** read_topol(const char* infile,
373 PreprocessingAtomTypes* atypes,
374 std::vector<MoleculeInformation>* molinfo,
375 std::unique_ptr<MoleculeInformation>* intermolecular_interactions,
376 gmx::ArrayRef<InteractionsOfType> interactions,
377 int* combination_rule,
381 std::vector<gmx_molblock_t>* molblock,
382 bool* ffParametrizedWithHBondConstraints,
385 bool usingFullRangeElectrostatics,
387 const gmx::MDLogger& logger)
391 char * pline = nullptr, **title = nullptr;
392 char line[STRLEN], errbuf[256], comb_str[256], nb_str[256];
394 char * dirstr, *dummy2;
395 int nrcopies, nscan, ncombs, ncopy;
396 double fLJ, fQQ, fPOW;
397 MoleculeInformation* mi0 = nullptr;
400 t_nbparam ** nbparam, **pair;
401 real fudgeLJ = -1; /* Multiplication factor to generate 1-4 from LJ */
402 bool bReadDefaults, bReadMolType, bGenPairs, bWarn_copy_A_B;
403 double qt = 0, qBt = 0; /* total charge */
404 int dcatt = -1, nmol_couple;
405 /* File handling variables */
409 char* tmp_line = nullptr;
410 char warn_buf[STRLEN];
411 const char* floating_point_arithmetic_tip =
412 "Total charge should normally be an integer. See\n"
413 "http://www.gromacs.org/Documentation/Floating_Point_Arithmetic\n"
414 "for discussion on how close it should be to an integer.\n";
415 /* We need to open the output file before opening the input file,
416 * because cpp_open_file can change the current working directory.
420 out = gmx_fio_fopen(outfile, "w");
427 /* open input file */
428 auto cpp_opts_return = cpp_opts(define, include, wi);
429 status = cpp_open_file(infile, &handle, cpp_opts_return);
432 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
435 /* some local variables */
436 DS_Init(&DS); /* directive stack */
437 d = Directive::d_invalid; /* first thing should be a directive */
438 nbparam = nullptr; /* The temporary non-bonded matrix */
439 pair = nullptr; /* The temporary pair interaction matrix */
440 std::vector<std::vector<gmx::ExclusionBlock>> exclusionBlocks;
443 *reppow = 12.0; /* Default value for repulsion power */
445 /* Init the number of CMAP torsion angles and grid spacing */
446 interactions[F_CMAP].cmakeGridSpacing = 0;
447 interactions[F_CMAP].cmapAngles = 0;
449 bWarn_copy_A_B = bFEP;
451 PreprocessingBondAtomType bondAtomType;
452 /* parse the actual file */
453 bReadDefaults = FALSE;
455 bReadMolType = FALSE;
460 status = cpp_read_line(&handle, STRLEN, line);
461 done = (status == eCPP_EOF);
464 if (status != eCPP_OK)
466 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
470 fprintf(out, "%s\n", line);
473 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
475 pline = gmx_strdup(line);
477 /* Strip trailing '\' from pline, if it exists */
479 if ((sl > 0) && (pline[sl - 1] == CONTINUE))
484 /* build one long line from several fragments - necessary for CMAP */
485 while (continuing(line))
487 status = cpp_read_line(&handle, STRLEN, line);
488 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
490 /* Since we depend on the '\' being present to continue to read, we copy line
491 * to a tmp string, strip the '\' from that string, and cat it to pline
493 tmp_line = gmx_strdup(line);
495 sl = strlen(tmp_line);
496 if ((sl > 0) && (tmp_line[sl - 1] == CONTINUE))
498 tmp_line[sl - 1] = ' ';
501 done = (status == eCPP_EOF);
504 if (status != eCPP_OK)
506 gmx_fatal(FARGS, "%s", cpp_error(&handle, status));
510 fprintf(out, "%s\n", line);
514 srenew(pline, strlen(pline) + strlen(tmp_line) + 1);
515 strcat(pline, tmp_line);
519 /* skip trailing and leading spaces and comment text */
520 strip_comment(pline);
523 /* if there is something left... */
524 if (static_cast<int>(strlen(pline)) > 0)
526 if (pline[0] == OPENDIR)
528 /* A directive on this line: copy the directive
529 * without the brackets into dirstr, then
530 * skip spaces and tabs on either side of directive
532 dirstr = gmx_strdup((pline + 1));
533 if ((dummy2 = strchr(dirstr, CLOSEDIR)) != nullptr)
539 if ((newd = str2dir(dirstr)) == Directive::d_invalid)
541 sprintf(errbuf, "Invalid directive %s", dirstr);
542 warning_error(wi, errbuf);
546 /* Directive found */
547 if (DS_Check_Order(DS, newd))
554 /* we should print here which directives should have
555 been present, and which actually are */
556 gmx_fatal(FARGS, "%s\nInvalid order for directive %s",
557 cpp_error(&handle, eCPP_SYNTAX), dir2str(newd));
558 /* d = Directive::d_invalid; */
561 if (d == Directive::d_intermolecular_interactions)
563 if (*intermolecular_interactions == nullptr)
565 /* We (mis)use the moleculetype processing
566 * to process the intermolecular interactions
567 * by making a "molecule" of the size of the system.
569 *intermolecular_interactions = std::make_unique<MoleculeInformation>();
570 mi0 = intermolecular_interactions->get();
572 make_atoms_sys(*molblock, *molinfo, &mi0->atoms);
578 else if (d != Directive::d_invalid)
580 /* Not a directive, just a plain string
581 * use a gigantic switch to decode,
582 * if there is a valid directive!
586 case Directive::d_defaults:
589 gmx_fatal(FARGS, "%s\nFound a second defaults directive.\n",
590 cpp_error(&handle, eCPP_SYNTAX));
592 bReadDefaults = TRUE;
593 nscan = sscanf(pline, "%s%s%s%lf%lf%lf", nb_str, comb_str, genpairs,
605 get_nbparm(nb_str, comb_str, &nb_funct, combination_rule, wi);
608 bGenPairs = (gmx::equalCaseInsensitive(genpairs, "Y", 1));
609 if (nb_funct != eNBF_LJ && bGenPairs)
612 "Generating pair parameters is only supported "
613 "with LJ non-bonded interactions");
629 nb_funct = ifunc_index(Directive::d_nonbond_params, nb_funct);
632 case Directive::d_atomtypes:
633 push_at(symtab, atypes, &bondAtomType, pline, nb_funct, &nbparam,
634 bGenPairs ? &pair : nullptr, wi);
637 case Directive::d_bondtypes:
638 push_bt(d, interactions, 2, nullptr, &bondAtomType, pline, wi);
640 case Directive::d_constrainttypes:
641 push_bt(d, interactions, 2, nullptr, &bondAtomType, pline, wi);
643 case Directive::d_pairtypes:
646 push_nbt(d, pair, atypes, pline, F_LJ14, wi);
650 push_bt(d, interactions, 2, atypes, nullptr, pline, wi);
653 case Directive::d_angletypes:
654 push_bt(d, interactions, 3, nullptr, &bondAtomType, pline, wi);
656 case Directive::d_dihedraltypes:
657 /* Special routine that can read both 2 and 4 atom dihedral definitions. */
658 push_dihedraltype(d, interactions, &bondAtomType, pline, wi);
661 case Directive::d_nonbond_params:
662 push_nbt(d, nbparam, atypes, pline, nb_funct, wi);
665 case Directive::d_implicit_genborn_params:
666 // Skip this line, so old topologies with
667 // GB parameters can be read.
670 case Directive::d_implicit_surface_params:
671 // Skip this line, so that any topologies
672 // with surface parameters can be read
673 // (even though these were never formally
677 case Directive::d_cmaptypes:
678 push_cmaptype(d, interactions, 5, atypes, &bondAtomType, pline, wi);
681 case Directive::d_moleculetype:
686 if (opts->couple_moltype != nullptr
687 && (opts->couple_lam0 == ecouplamNONE || opts->couple_lam0 == ecouplamQ
688 || opts->couple_lam1 == ecouplamNONE
689 || opts->couple_lam1 == ecouplamQ))
691 dcatt = add_atomtype_decoupled(symtab, atypes, &nbparam,
692 bGenPairs ? &pair : nullptr);
694 ntype = atypes->size();
695 ncombs = (ntype * (ntype + 1)) / 2;
696 generate_nbparams(*combination_rule, nb_funct,
697 &(interactions[nb_funct]), atypes, wi);
698 ncopy = copy_nbparams(nbparam, nb_funct, &(interactions[nb_funct]), ntype);
701 .appendTextFormatted(
702 "Generated %d of the %d non-bonded parameter "
704 ncombs - ncopy, ncombs);
705 free_nbparam(nbparam, ntype);
708 gen_pairs((interactions[nb_funct]), &(interactions[F_LJ14]),
709 fudgeLJ, *combination_rule);
710 ncopy = copy_nbparams(pair, nb_funct, &(interactions[F_LJ14]), ntype);
713 .appendTextFormatted(
714 "Generated %d of the %d 1-4 parameter "
716 ncombs - ncopy, ncombs);
717 free_nbparam(pair, ntype);
719 /* Copy GBSA parameters to atomtype array? */
724 push_molt(symtab, molinfo, pline, wi);
725 exclusionBlocks.emplace_back();
726 mi0 = &molinfo->back();
727 mi0->atoms.haveMass = TRUE;
728 mi0->atoms.haveCharge = TRUE;
729 mi0->atoms.haveType = TRUE;
730 mi0->atoms.haveBState = TRUE;
731 mi0->atoms.havePdbInfo = FALSE;
734 case Directive::d_atoms:
735 push_atom(symtab, &(mi0->atoms), atypes, pline, wi);
738 case Directive::d_pairs:
741 "Need to have a valid MoleculeInformation object to work on");
742 push_bond(d, interactions, mi0->interactions, &(mi0->atoms), atypes,
743 pline, FALSE, bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
745 case Directive::d_pairs_nb:
748 "Need to have a valid MoleculeInformation object to work on");
749 push_bond(d, interactions, mi0->interactions, &(mi0->atoms), atypes,
750 pline, FALSE, FALSE, 1.0, bZero, &bWarn_copy_A_B, wi);
753 case Directive::d_vsites2:
754 case Directive::d_vsites3:
755 case Directive::d_vsites4:
756 case Directive::d_bonds:
757 case Directive::d_angles:
758 case Directive::d_constraints:
759 case Directive::d_settles:
760 case Directive::d_position_restraints:
761 case Directive::d_angle_restraints:
762 case Directive::d_angle_restraints_z:
763 case Directive::d_distance_restraints:
764 case Directive::d_orientation_restraints:
765 case Directive::d_dihedral_restraints:
766 case Directive::d_dihedrals:
767 case Directive::d_polarization:
768 case Directive::d_water_polarization:
769 case Directive::d_thole_polarization:
772 "Need to have a valid MoleculeInformation object to work on");
773 push_bond(d, interactions, mi0->interactions, &(mi0->atoms), atypes,
774 pline, TRUE, bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
776 case Directive::d_cmap:
779 "Need to have a valid MoleculeInformation object to work on");
780 push_cmap(d, interactions, mi0->interactions, &(mi0->atoms), atypes, pline, wi);
783 case Directive::d_vsitesn:
786 "Need to have a valid MoleculeInformation object to work on");
787 push_vsitesn(d, mi0->interactions, &(mi0->atoms), pline, wi);
789 case Directive::d_exclusions:
790 GMX_ASSERT(!exclusionBlocks.empty(),
791 "exclusionBlocks must always be allocated so exclusions can "
793 if (exclusionBlocks.back().empty())
795 GMX_RELEASE_ASSERT(mi0,
796 "Need to have a valid MoleculeInformation "
797 "object to work on");
798 exclusionBlocks.back().resize(mi0->atoms.nr);
800 push_excl(pline, exclusionBlocks.back(), wi);
802 case Directive::d_system:
804 title = put_symtab(symtab, pline);
806 case Directive::d_molecules:
811 push_mol(*molinfo, pline, &whichmol, &nrcopies, wi);
812 mi0 = &((*molinfo)[whichmol]);
813 molblock->resize(molblock->size() + 1);
814 molblock->back().type = whichmol;
815 molblock->back().nmol = nrcopies;
817 bCouple = (opts->couple_moltype != nullptr
818 && (gmx_strcasecmp("system", opts->couple_moltype) == 0
819 || strcmp(*(mi0->name), opts->couple_moltype) == 0));
822 nmol_couple += nrcopies;
825 if (mi0->atoms.nr == 0)
827 gmx_fatal(FARGS, "Molecule type '%s' contains no atoms", *mi0->name);
831 .appendTextFormatted(
832 "Excluding %d bonded neighbours molecule type '%s'",
833 mi0->nrexcl, *mi0->name);
834 sum_q(&mi0->atoms, nrcopies, &qt, &qBt);
835 if (!mi0->bProcessed)
837 generate_excl(mi0->nrexcl, mi0->atoms.nr, mi0->interactions, &(mi0->excls));
838 gmx::mergeExclusions(&(mi0->excls), exclusionBlocks[whichmol]);
839 make_shake(mi0->interactions, &mi0->atoms, opts->nshake, logger);
843 convert_moltype_couple(mi0, dcatt, *fudgeQQ, opts->couple_lam0,
844 opts->couple_lam1, opts->bCoupleIntra,
845 nb_funct, &(interactions[nb_funct]), wi);
847 stupid_fill_block(&mi0->mols, mi0->atoms.nr, TRUE);
848 mi0->bProcessed = TRUE;
853 GMX_LOG(logger.warning)
855 .appendTextFormatted("case: %d", static_cast<int>(d));
856 gmx_incons("unknown directive");
865 // Check that all strings defined with -D were used when processing topology
866 std::string unusedDefineWarning = checkAndWarnForUnusedDefines(*handle);
867 if (!unusedDefineWarning.empty())
869 warning(wi, unusedDefineWarning);
872 sfree(cpp_opts_return);
879 /* List of GROMACS define names for force fields that have been
880 * parametrized using constraints involving hydrogens only.
882 * We should avoid hardcoded names, but this is hopefully only
883 * needed temparorily for discouraging use of constraints=all-bonds.
885 const std::array<std::string, 3> ffDefines = { "_FF_AMBER", "_FF_CHARMM", "_FF_OPLSAA" };
886 *ffParametrizedWithHBondConstraints = false;
887 for (const std::string& ffDefine : ffDefines)
889 if (cpp_find_define(&handle, ffDefine))
891 *ffParametrizedWithHBondConstraints = true;
895 if (cpp_find_define(&handle, "_FF_GROMOS96") != nullptr)
898 "The GROMOS force fields have been parametrized with a physically incorrect "
899 "multiple-time-stepping scheme for a twin-range cut-off. When used with "
900 "a single-range cut-off (or a correct Trotter multiple-time-stepping scheme), "
901 "physical properties, such as the density, might differ from the intended values. "
902 "Since there are researchers actively working on validating GROMOS with modern "
903 "integrators we have not yet removed the GROMOS force fields, but you should be "
904 "aware of these issues and check if molecules in your system are affected before "
906 "Further information is available at https://redmine.gromacs.org/issues/2884 , "
907 "and a longer explanation of our decision to remove physically incorrect "
909 "can be found at https://doi.org/10.26434/chemrxiv.11474583.v1 .");
914 if (opts->couple_moltype)
916 if (nmol_couple == 0)
918 gmx_fatal(FARGS, "Did not find any molecules of type '%s' for coupling", opts->couple_moltype);
922 .appendTextFormatted("Coupling %d copies of molecule type '%s'", nmol_couple,
923 opts->couple_moltype);
926 /* this is not very clean, but fixes core dump on empty system name */
929 title = put_symtab(symtab, "");
934 sprintf(warn_buf, "System has non-zero total charge: %.6f\n%s\n", qt, floating_point_arithmetic_tip);
935 warning_note(wi, warn_buf);
937 if (fabs(qBt) > 1e-4 && !gmx_within_tol(qBt, qt, 1e-6))
939 sprintf(warn_buf, "State B has non-zero total charge: %.6f\n%s\n", qBt,
940 floating_point_arithmetic_tip);
941 warning_note(wi, warn_buf);
943 if (usingFullRangeElectrostatics && (fabs(qt) > 1e-4 || fabs(qBt) > 1e-4))
946 "You are using Ewald electrostatics in a system with net charge. This can lead to "
947 "severe artifacts, such as ions moving into regions with low dielectric, due to "
948 "the uniform background charge. We suggest to neutralize your system with counter "
949 "ions, possibly in combination with a physiological salt concentration.");
950 please_cite(stdout, "Hub2014a");
955 if (*intermolecular_interactions != nullptr)
957 sfree(intermolecular_interactions->get()->atoms.atom);
963 char** do_top(bool bVerbose,
965 const char* topppfile,
969 gmx::ArrayRef<InteractionsOfType> interactions,
970 int* combination_rule,
971 double* repulsion_power,
973 PreprocessingAtomTypes* atypes,
974 std::vector<MoleculeInformation>* molinfo,
975 std::unique_ptr<MoleculeInformation>* intermolecular_interactions,
976 const t_inputrec* ir,
977 std::vector<gmx_molblock_t>* molblock,
978 bool* ffParametrizedWithHBondConstraints,
980 const gmx::MDLogger& logger)
982 /* Tmpfile might contain a long path */
997 GMX_LOG(logger.info).asParagraph().appendTextFormatted("processing topology...");
999 title = read_topol(topfile, tmpfile, opts->define, opts->include, symtab, atypes, molinfo,
1000 intermolecular_interactions, interactions, combination_rule, repulsion_power,
1001 opts, fudgeQQ, molblock, ffParametrizedWithHBondConstraints,
1002 ir->efep != efepNO, bZero, EEL_FULL(ir->coulombtype), wi, logger);
1004 if ((*combination_rule != eCOMB_GEOMETRIC) && (ir->vdwtype == evdwUSER))
1007 "Using sigma/epsilon based combination rules with"
1008 " user supplied potential function may produce unwanted"
1016 * Generate exclusion lists for QM/MM.
1018 * This routine updates the exclusion lists for QM atoms in order to include all other QM
1019 * atoms of this molecule. Moreover, this routine replaces bonds between QM atoms with
1020 * CONNBOND and, when MiMiC is not used, removes bonded interactions between QM and link atoms.
1021 * Finally, in case if MiMiC QM/MM is used - charges of QM atoms are set to 0
1023 * \param[in,out] molt molecule type with QM atoms
1024 * \param[in] grpnr group informatio
1025 * \param[in,out] ir input record
1026 * \param[in,out] qmmmMode QM/MM mode switch: original/MiMiC
1027 * \param[in] logger Handle to logging interface.
1029 static void generate_qmexcl_moltype(gmx_moltype_t* molt,
1030 const unsigned char* grpnr,
1032 GmxQmmmMode qmmmMode,
1033 const gmx::MDLogger& logger)
1035 /* This routine expects molt->ilist to be of size F_NRE and ordered. */
1037 /* generates the exclusions between the individual QM atoms, as
1038 * these interactions should be handled by the QM subroutines and
1039 * not by the gromacs routines
1041 int qm_max = 0, qm_nr = 0, link_nr = 0, link_max = 0;
1042 int * qm_arr = nullptr, *link_arr = nullptr;
1043 bool *bQMMM, *blink;
1045 /* First we search and select the QM atoms in an qm_arr array that
1046 * we use to create the exclusions.
1048 * we take the possibility into account that a user has defined more
1049 * than one QM group:
1051 * for that we also need to do this an ugly work-about just in case
1052 * the QM group contains the entire system...
1055 /* we first search for all the QM atoms and put them in an array
1057 for (int j = 0; j < ir->opts.ngQM; j++)
1059 for (int i = 0; i < molt->atoms.nr; i++)
1061 if (qm_nr >= qm_max)
1064 srenew(qm_arr, qm_max);
1066 if ((grpnr ? grpnr[i] : 0) == j)
1068 qm_arr[qm_nr++] = i;
1069 molt->atoms.atom[i].q = 0.0;
1070 molt->atoms.atom[i].qB = 0.0;
1074 /* bQMMM[..] is an array containin TRUE/FALSE for atoms that are
1075 * QM/not QM. We first set all elements to false. Afterwards we use
1076 * the qm_arr to change the elements corresponding to the QM atoms
1079 snew(bQMMM, molt->atoms.nr);
1080 for (int i = 0; i < molt->atoms.nr; i++)
1084 for (int i = 0; i < qm_nr; i++)
1086 bQMMM[qm_arr[i]] = TRUE;
1089 /* We remove all bonded interactions (i.e. bonds,
1090 * angles, dihedrals, 1-4's), involving the QM atoms. The way they
1091 * are removed is as follows: if the interaction invloves 2 atoms,
1092 * it is removed if both atoms are QMatoms. If it involves 3 atoms,
1093 * it is removed if at least two of the atoms are QM atoms, if the
1094 * interaction involves 4 atoms, it is removed if there are at least
1095 * 2 QM atoms. Since this routine is called once before any forces
1096 * are computed, the top->idef.il[N].iatom[] array (see idef.h) can
1097 * be rewritten at this poitn without any problem. 25-9-2002 */
1099 /* first check whether we already have CONNBONDS.
1100 * Note that if we don't, we don't add a param entry and set ftype=0,
1101 * which is ok, since CONNBONDS does not use parameters.
1103 int ftype_connbond = 0;
1104 int ind_connbond = 0;
1105 if (molt->ilist[F_CONNBONDS].size() != 0)
1107 GMX_LOG(logger.info)
1109 .appendTextFormatted("nr. of CONNBONDS present already: %d",
1110 molt->ilist[F_CONNBONDS].size() / 3);
1111 ftype_connbond = molt->ilist[F_CONNBONDS].iatoms[0];
1112 ind_connbond = molt->ilist[F_CONNBONDS].size();
1114 /* now we delete all bonded interactions, except the ones describing
1115 * a chemical bond. These are converted to CONNBONDS
1117 for (int ftype = 0; ftype < F_NRE; ftype++)
1119 if (!(interaction_function[ftype].flags & IF_BOND) || ftype == F_CONNBONDS)
1123 int nratoms = interaction_function[ftype].nratoms;
1125 while (j < molt->ilist[ftype].size())
1131 /* Remove an interaction between two atoms when both are
1132 * in the QM region. Note that we don't have to worry about
1133 * link atoms here, as they won't have 2-atom interactions.
1135 int a1 = molt->ilist[ftype].iatoms[1 + j + 0];
1136 int a2 = molt->ilist[ftype].iatoms[1 + j + 1];
1137 bexcl = (bQMMM[a1] && bQMMM[a2]);
1138 /* A chemical bond between two QM atoms will be copied to
1139 * the F_CONNBONDS list, for reasons mentioned above.
1141 if (bexcl && IS_CHEMBOND(ftype))
1143 InteractionList& ilist = molt->ilist[F_CONNBONDS];
1144 ilist.iatoms.resize(ind_connbond + 3);
1145 ilist.iatoms[ind_connbond++] = ftype_connbond;
1146 ilist.iatoms[ind_connbond++] = a1;
1147 ilist.iatoms[ind_connbond++] = a2;
1152 /* MM interactions have to be excluded if they are included
1153 * in the QM already. Because we use a link atom (H atom)
1154 * when the QM/MM boundary runs through a chemical bond, this
1155 * means that as long as one atom is MM, we still exclude,
1156 * as the interaction is included in the QM via:
1157 * QMatom1-QMatom2-QMatom-3-Linkatom.
1160 for (int jj = j + 1; jj < j + 1 + nratoms; jj++)
1162 if (bQMMM[molt->ilist[ftype].iatoms[jj]])
1168 /* MiMiC treats link atoms as quantum atoms - therefore
1169 * we do not need do additional exclusions here */
1170 if (qmmmMode == GmxQmmmMode::GMX_QMMM_MIMIC)
1172 bexcl = numQmAtoms == nratoms;
1176 bexcl = (numQmAtoms >= nratoms - 1);
1179 if (bexcl && ftype == F_SETTLE)
1182 "Can not apply QM to molecules with SETTLE, replace the moleculetype "
1183 "using QM and SETTLE by one without SETTLE");
1188 /* since the interaction involves QM atoms, these should be
1189 * removed from the MM ilist
1191 InteractionList& ilist = molt->ilist[ftype];
1192 for (int k = j; k < ilist.size() - (nratoms + 1); k++)
1194 ilist.iatoms[k] = ilist.iatoms[k + (nratoms + 1)];
1196 ilist.iatoms.resize(ilist.size() - (nratoms + 1));
1200 j += nratoms + 1; /* the +1 is for the functype */
1204 /* Now, we search for atoms bonded to a QM atom because we also want
1205 * to exclude their nonbonded interactions with the QM atoms. The
1206 * reason for this is that this interaction is accounted for in the
1207 * linkatoms interaction with the QMatoms and would be counted
1210 if (qmmmMode != GmxQmmmMode::GMX_QMMM_MIMIC)
1212 for (int i = 0; i < F_NRE; i++)
1217 while (j < molt->ilist[i].size())
1219 int a1 = molt->ilist[i].iatoms[j + 1];
1220 int a2 = molt->ilist[i].iatoms[j + 2];
1221 if ((bQMMM[a1] && !bQMMM[a2]) || (!bQMMM[a1] && bQMMM[a2]))
1223 if (link_nr >= link_max)
1226 srenew(link_arr, link_max);
1230 link_arr[link_nr++] = a2;
1234 link_arr[link_nr++] = a1;
1242 snew(blink, molt->atoms.nr);
1243 for (int i = 0; i < molt->atoms.nr; i++)
1248 if (qmmmMode != GmxQmmmMode::GMX_QMMM_MIMIC)
1250 for (int i = 0; i < link_nr; i++)
1252 blink[link_arr[i]] = TRUE;
1255 /* creating the exclusion block for the QM atoms. Each QM atom has
1256 * as excluded elements all the other QMatoms (and itself).
1259 qmexcl.nr = molt->atoms.nr;
1260 qmexcl.nra = qm_nr * (qm_nr + link_nr) + link_nr * qm_nr;
1261 snew(qmexcl.index, qmexcl.nr + 1);
1262 snew(qmexcl.a, qmexcl.nra);
1264 for (int i = 0; i < qmexcl.nr; i++)
1266 qmexcl.index[i] = j;
1269 for (int k = 0; k < qm_nr; k++)
1271 qmexcl.a[k + j] = qm_arr[k];
1273 for (int k = 0; k < link_nr; k++)
1275 qmexcl.a[qm_nr + k + j] = link_arr[k];
1277 j += (qm_nr + link_nr);
1281 for (int k = 0; k < qm_nr; k++)
1283 qmexcl.a[k + j] = qm_arr[k];
1288 qmexcl.index[qmexcl.nr] = j;
1290 /* and merging with the exclusions already present in sys.
1293 std::vector<gmx::ExclusionBlock> qmexcl2(molt->atoms.nr);
1294 gmx::blockaToExclusionBlocks(&qmexcl, qmexcl2);
1295 gmx::mergeExclusions(&(molt->excls), qmexcl2);
1297 /* Finally, we also need to get rid of the pair interactions of the
1298 * classical atom bonded to the boundary QM atoms with the QMatoms,
1299 * as this interaction is already accounted for by the QM, so also
1300 * here we run the risk of double counting! We proceed in a similar
1301 * way as we did above for the other bonded interactions: */
1302 for (int i = F_LJ14; i < F_COUL14; i++)
1304 int nratoms = interaction_function[i].nratoms;
1306 while (j < molt->ilist[i].size())
1308 int a1 = molt->ilist[i].iatoms[j + 1];
1309 int a2 = molt->ilist[i].iatoms[j + 2];
1311 ((bQMMM[a1] && bQMMM[a2]) || (blink[a1] && bQMMM[a2]) || (bQMMM[a1] && blink[a2]));
1314 /* since the interaction involves QM atoms, these should be
1315 * removed from the MM ilist
1317 InteractionList& ilist = molt->ilist[i];
1318 for (int k = j; k < ilist.size() - (nratoms + 1); k++)
1320 ilist.iatoms[k] = ilist.iatoms[k + (nratoms + 1)];
1322 ilist.iatoms.resize(ilist.size() - (nratoms + 1));
1326 j += nratoms + 1; /* the +1 is for the functype */
1335 } /* generate_qmexcl */
1337 void generate_qmexcl(gmx_mtop_t* sys, t_inputrec* ir, warninp* wi, GmxQmmmMode qmmmMode, const gmx::MDLogger& logger)
1339 /* This routine expects molt->molt[m].ilist to be of size F_NRE and ordered.
1342 unsigned char* grpnr;
1343 int mol, nat_mol, nr_mol_with_qm_atoms = 0;
1344 gmx_molblock_t* molb;
1346 int index_offset = 0;
1349 grpnr = sys->groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].data();
1351 for (size_t mb = 0; mb < sys->molblock.size(); mb++)
1353 molb = &sys->molblock[mb];
1354 nat_mol = sys->moltype[molb->type].atoms.nr;
1355 for (mol = 0; mol < molb->nmol; mol++)
1358 for (int i = 0; i < nat_mol; i++)
1360 if ((grpnr ? grpnr[i] : 0) < (ir->opts.ngQM))
1369 nr_mol_with_qm_atoms++;
1372 /* We need to split this molblock */
1375 /* Split the molblock at this molecule */
1376 auto pos = sys->molblock.begin() + mb + 1;
1377 sys->molblock.insert(pos, sys->molblock[mb]);
1378 sys->molblock[mb].nmol = mol;
1379 sys->molblock[mb + 1].nmol -= mol;
1381 molb = &sys->molblock[mb];
1385 /* Split the molblock after this molecule */
1386 auto pos = sys->molblock.begin() + mb + 1;
1387 sys->molblock.insert(pos, sys->molblock[mb]);
1388 molb = &sys->molblock[mb];
1389 sys->molblock[mb].nmol = 1;
1390 sys->molblock[mb + 1].nmol -= 1;
1393 /* Create a copy of a moltype for a molecule
1394 * containing QM atoms and append it in the end of the list
1396 std::vector<gmx_moltype_t> temp(sys->moltype.size());
1397 for (size_t i = 0; i < sys->moltype.size(); ++i)
1399 copy_moltype(&sys->moltype[i], &temp[i]);
1401 sys->moltype.resize(sys->moltype.size() + 1);
1402 for (size_t i = 0; i < temp.size(); ++i)
1404 copy_moltype(&temp[i], &sys->moltype[i]);
1406 copy_moltype(&sys->moltype[molb->type], &sys->moltype.back());
1407 /* Copy the exclusions to a new array, since this is the only
1408 * thing that needs to be modified for QMMM.
1410 sys->moltype.back().excls = sys->moltype[molb->type].excls;
1411 /* Set the molecule type for the QMMM molblock */
1412 molb->type = sys->moltype.size() - 1;
1414 generate_qmexcl_moltype(&sys->moltype[molb->type], grpnr, ir, qmmmMode, logger);
1420 index_offset += nat_mol;
1423 if (qmmmMode == GmxQmmmMode::GMX_QMMM_ORIGINAL && nr_mol_with_qm_atoms > 1)
1425 /* generate a warning is there are QM atoms in different
1426 * topologies. In this case it is not possible at this stage to
1427 * mutualy exclude the non-bonded interactions via the
1428 * exclusions (AFAIK). Instead, the user is advised to use the
1429 * energy group exclusions in the mdp file
1432 "\nThe QM subsystem is divided over multiple topologies. "
1433 "The mutual non-bonded interactions cannot be excluded. "
1434 "There are two ways to achieve this:\n\n"
1435 "1) merge the topologies, such that the atoms of the QM "
1436 "subsystem are all present in one single topology file. "
1437 "In this case this warning will dissappear\n\n"
1438 "2) exclude the non-bonded interactions explicitly via the "
1439 "energygrp-excl option in the mdp file. if this is the case "
1440 "this warning may be ignored"