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, 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.
42 #include <sys/types.h>
49 #include "gromacs/fileio/futil.h"
54 #include "gromacs/fileio/gmxfio.h"
61 #include "gmx_fatal.h"
63 #include "vsite_parm.h"
65 #include "grompp-impl.h"
69 #include "gpp_nextnb.h"
73 #include "gpp_bond_atomtype.h"
75 #include "gromacs/math/utilities.h"
77 #define CPPMARK '#' /* mark from cpp */
78 #define OPENDIR '[' /* starting sign for directive */
79 #define CLOSEDIR ']' /* ending sign for directive */
82 printf("line: %d, maxavail: %d\n", __LINE__, maxavail()); \
85 static void free_nbparam(t_nbparam **param, int nr)
89 for (i = 0; i < nr; i++)
96 static int copy_nbparams(t_nbparam **param, int ftype, t_params *plist, int nr)
104 for (i = 0; i < nr; i++)
106 for (j = 0; j <= i; j++)
108 if (param[i][j].bSet)
110 for (f = 0; f < nrfp; f++)
112 plist->param[nr*i+j].c[f] = param[i][j].c[f];
113 plist->param[nr*j+i].c[f] = param[i][j].c[f];
123 static void gen_pairs(t_params *nbs, t_params *pairs, real fudge, int comb)
125 int i, j, ntp, nrfp, nrfpA, nrfpB, nnn;
130 nrfpA = interaction_function[F_LJ14].nrfpA;
131 nrfpB = interaction_function[F_LJ14].nrfpB;
134 if ((nrfp != nrfpA) || (nrfpA != nrfpB))
136 gmx_incons("Number of force parameters in gen_pairs wrong");
139 fprintf(stderr, "Generating 1-4 interactions: fudge = %g\n", fudge);
142 fprintf(debug, "Fudge factor for 1-4 interactions: %g\n", fudge);
143 fprintf(debug, "Holy Cow! there are %d types\n", ntp);
145 snew(pairs->param, pairs->nr);
146 for (i = 0; (i < ntp); i++)
149 pairs->param[i].a[0] = i / nnn;
150 pairs->param[i].a[1] = i % nnn;
151 /* Copy normal and FEP parameters and multiply by fudge factor */
155 for (j = 0; (j < nrfp); j++)
157 /* If we are using sigma/epsilon values, only the epsilon values
158 * should be scaled, but not sigma.
159 * The sigma values have even indices 0,2, etc.
161 if ((comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) && (j%2 == 0))
170 pairs->param[i].c[j] = scaling*nbs->param[i].c[j];
171 pairs->param[i].c[nrfp+j] = scaling*nbs->param[i].c[j];
176 double check_mol(gmx_mtop_t *mtop, warninp_t wi)
179 int i, mb, nmol, ri, pt;
184 /* Check mass and charge */
187 for (mb = 0; mb < mtop->nmoltype; mb++)
189 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
190 nmol = mtop->molblock[mb].nmol;
191 for (i = 0; (i < atoms->nr); i++)
193 q += nmol*atoms->atom[i].q;
194 m = atoms->atom[i].m;
195 pt = atoms->atom[i].ptype;
196 /* If the particle is an atom or a nucleus it must have a mass,
197 * else, if it is a shell, a vsite or a bondshell it can have mass zero
199 if ((m <= 0.0) && ((pt == eptAtom) || (pt == eptNucleus)))
201 ri = atoms->atom[i].resind;
202 sprintf(buf, "atom %s (Res %s-%d) has mass %g\n",
203 *(atoms->atomname[i]),
204 *(atoms->resinfo[ri].name),
205 atoms->resinfo[ri].nr,
207 warning_error(wi, buf);
210 if ((m != 0) && (pt == eptVSite))
212 ri = atoms->atom[i].resind;
213 sprintf(buf, "virtual site %s (Res %s-%d) has non-zero mass %g\n"
214 " Check your topology.\n",
215 *(atoms->atomname[i]),
216 *(atoms->resinfo[ri].name),
217 atoms->resinfo[ri].nr,
219 warning_error(wi, buf);
220 /* The following statements make LINCS break! */
221 /* atoms->atom[i].m=0; */
228 static void sum_q(t_atoms *atoms, int n, double *qt, double *qBt)
236 for (i = 0; i < atoms->nr; i++)
238 qmolA += atoms->atom[i].q;
239 qmolB += atoms->atom[i].qB;
241 /* Unfortunately an absolute comparison,
242 * but this avoids unnecessary warnings and gmx-users mails.
244 if (fabs(qmolA) >= 1e-6 || fabs(qmolB) >= 1e-6)
251 static void get_nbparm(char *nb_str, char *comb_str, int *nb, int *comb,
255 char warn_buf[STRLEN];
258 for (i = 1; (i < eNBF_NR); i++)
260 if (gmx_strcasecmp(nb_str, enbf_names[i]) == 0)
267 *nb = strtol(nb_str, NULL, 10);
269 if ((*nb < 1) || (*nb >= eNBF_NR))
271 sprintf(warn_buf, "Invalid nonbond function selector '%s' using %s",
272 nb_str, enbf_names[1]);
273 warning_error(wi, warn_buf);
277 for (i = 1; (i < eCOMB_NR); i++)
279 if (gmx_strcasecmp(comb_str, ecomb_names[i]) == 0)
286 *comb = strtol(comb_str, NULL, 10);
288 if ((*comb < 1) || (*comb >= eCOMB_NR))
290 sprintf(warn_buf, "Invalid combination rule selector '%s' using %s",
291 comb_str, ecomb_names[1]);
292 warning_error(wi, warn_buf);
297 static char ** cpp_opts(const char *define, const char *include,
302 const char *cppadds[2];
303 char **cppopts = NULL;
304 const char *option[2] = { "-D", "-I" };
305 const char *nopt[2] = { "define", "include" };
309 char warn_buf[STRLEN];
312 cppadds[1] = include;
313 for (n = 0; (n < 2); n++)
320 while ((*ptr != '\0') && isspace(*ptr))
325 while ((*rptr != '\0') && !isspace(*rptr))
333 strncpy(buf, ptr, len);
334 if (strstr(ptr, option[n]) != ptr)
336 set_warning_line(wi, "mdp file", -1);
337 sprintf(warn_buf, "Malformed %s option %s", nopt[n], buf);
338 warning(wi, warn_buf);
342 srenew(cppopts, ++ncppopts);
343 cppopts[ncppopts-1] = strdup(buf);
351 srenew(cppopts, ++ncppopts);
352 cppopts[ncppopts-1] = NULL;
359 find_gb_bondlength(t_params *plist, int ai, int aj, real *length)
366 for (i = 0; i < F_NRE && !found; i++)
370 for (j = 0; j < plist[i].nr; j++)
372 a1 = plist[i].param[j].a[0];
373 a2 = plist[i].param[j].a[1];
375 if ( (a1 == ai && a2 == aj) || (a1 == aj && a2 == ai))
377 /* Equilibrium bond distance */
378 *length = plist[i].param[j].c[0];
391 find_gb_anglelength(t_params *plist, int ai, int ak, real *length)
393 int i, j, a1, a2, a3;
396 int status, status1, status2;
400 for (i = 0; i < F_NRE && !found; i++)
404 for (j = 0; j < plist[i].nr; j++)
406 a1 = plist[i].param[j].a[0];
407 a2 = plist[i].param[j].a[1];
408 a3 = plist[i].param[j].a[2];
410 /* We dont care what the middle atom is, but use it below */
411 if ( (a1 == ai && a3 == ak) || (a1 == ak && a3 == ai) )
413 /* Equilibrium bond distance */
414 a123 = plist[i].param[j].c[0];
415 /* Use middle atom to find reference distances r12 and r23 */
416 status1 = find_gb_bondlength(plist, a1, a2, &r12);
417 status2 = find_gb_bondlength(plist, a2, a3, &r23);
419 if (status1 == 0 && status2 == 0)
421 /* cosine theorem to get r13 */
422 *length = sqrt(r12*r12+r23*r23-(2*r12*r23*cos(a123/RAD2DEG)));
435 generate_gb_exclusion_interactions(t_molinfo *mi, gpp_atomtype_t atype, t_nextnb *nnb)
437 int i, j, k, n, ai, aj, ti, tj;
443 real radiusi, radiusj;
444 real gb_radiusi, gb_radiusj;
445 real param_c2, param_c4;
451 for (n = 1; n <= nnb->nrex; n++)
466 /* Put all higher-order exclusions into 1,4 list so we dont miss them */
473 for (ai = 0; ai < nnb->nr; ai++)
475 ti = at->atom[ai].type;
476 radiusi = get_atomtype_radius(ti, atype);
477 gb_radiusi = get_atomtype_gb_radius(ti, atype);
479 for (j = 0; j < nnb->nrexcl[ai][n]; j++)
481 aj = nnb->a[ai][n][j];
483 /* Only add the interactions once */
486 tj = at->atom[aj].type;
487 radiusj = get_atomtype_radius(tj, atype);
488 gb_radiusj = get_atomtype_gb_radius(tj, atype);
490 /* There is an exclusion of type "ftype" between atoms ai and aj */
494 /* Reference distance, not used for 1-4 interactions */
498 if (find_gb_bondlength(plist, ai, aj, &distance) != 0)
500 gmx_fatal(FARGS, "Cannot find bond length for atoms %d-%d", ai, aj);
504 if (find_gb_anglelength(plist, ai, aj, &distance) != 0)
506 gmx_fatal(FARGS, "Cannot find length for atoms %d-%d involved in angle", ai, aj);
513 /* Assign GB parameters */
515 param.c[0] = radiusi+radiusj;
516 /* Reference distance distance */
517 param.c[1] = distance;
518 /* Still parameter */
519 param.c[2] = param_c2;
521 param.c[3] = gb_radiusi+gb_radiusj;
523 param.c[4] = param_c4;
525 /* Add it to the parameter list */
526 add_param_to_list(&plist[ftype], ¶m);
537 static char **read_topol(const char *infile, const char *outfile,
538 const char *define, const char *include,
540 gpp_atomtype_t atype,
544 int *combination_rule,
549 gmx_molblock_t **molblock,
556 int i, sl, nb_funct, comb;
557 char *pline = NULL, **title = NULL;
558 char line[STRLEN], errbuf[256], comb_str[256], nb_str[256];
560 char *dirstr, *dummy2;
561 int nrcopies, nmol, nmolb = 0, nscan, ncombs, ncopy;
562 double fLJ, fQQ, fPOW;
563 gmx_molblock_t *molb = NULL;
564 t_topology *block = NULL;
565 t_molinfo *mi0 = NULL;
568 t_nbparam **nbparam, **pair;
570 real fudgeLJ = -1; /* Multiplication factor to generate 1-4 from LJ */
571 gmx_bool bReadDefaults, bReadMolType, bGenPairs, bWarn_copy_A_B;
572 double qt = 0, qBt = 0; /* total charge */
573 t_bond_atomtype batype;
575 int dcatt = -1, nmol_couple;
576 /* File handling variables */
579 char *tmp_line = NULL;
580 char warn_buf[STRLEN];
581 const char *floating_point_arithmetic_tip =
582 "Total charge should normally be an integer. See\n"
583 "http://www.gromacs.org/Documentation/Floating_Point_Arithmetic\n"
584 "for discussion on how close it should be to an integer.\n";
585 /* We need to open the output file before opening the input file,
586 * because cpp_open_file can change the current working directory.
590 out = gmx_fio_fopen(outfile, "w");
597 /* open input file */
598 status = cpp_open_file(infile, &handle, cpp_opts(define, include, wi));
601 gmx_fatal(FARGS, cpp_error(&handle, status));
604 /* some local variables */
605 DS_Init(&DS); /* directive stack */
606 nmol = 0; /* no molecules yet... */
607 d = d_invalid; /* first thing should be a directive */
608 nbparam = NULL; /* The temporary non-bonded matrix */
609 pair = NULL; /* The temporary pair interaction matrix */
610 block2 = NULL; /* the extra exclusions */
612 *reppow = 12.0; /* Default value for repulsion power */
616 /* Init the number of CMAP torsion angles and grid spacing */
617 plist->grid_spacing = 0;
620 bWarn_copy_A_B = bFEP;
622 batype = init_bond_atomtype();
623 /* parse the actual file */
624 bReadDefaults = FALSE;
626 bReadMolType = FALSE;
631 status = cpp_read_line(&handle, STRLEN, line);
632 done = (status == eCPP_EOF);
635 if (status != eCPP_OK)
637 gmx_fatal(FARGS, cpp_error(&handle, status));
641 fprintf(out, "%s\n", line);
644 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
646 pline = strdup(line);
648 /* Strip trailing '\' from pline, if it exists */
650 if ((sl > 0) && (pline[sl-1] == CONTINUE))
655 /* build one long line from several fragments - necessary for CMAP */
656 while (continuing(line))
658 status = cpp_read_line(&handle, STRLEN, line);
659 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
661 /* Since we depend on the '\' being present to continue to read, we copy line
662 * to a tmp string, strip the '\' from that string, and cat it to pline
664 tmp_line = strdup(line);
666 sl = strlen(tmp_line);
667 if ((sl > 0) && (tmp_line[sl-1] == CONTINUE))
669 tmp_line[sl-1] = ' ';
672 done = (status == eCPP_EOF);
675 if (status != eCPP_OK)
677 gmx_fatal(FARGS, cpp_error(&handle, status));
681 fprintf(out, "%s\n", line);
685 srenew(pline, strlen(pline)+strlen(tmp_line)+1);
686 strcat(pline, tmp_line);
690 /* skip trailing and leading spaces and comment text */
691 strip_comment (pline);
694 /* if there is something left... */
695 if ((int)strlen(pline) > 0)
697 if (pline[0] == OPENDIR)
699 /* A directive on this line: copy the directive
700 * without the brackets into dirstr, then
701 * skip spaces and tabs on either side of directive
703 dirstr = strdup((pline+1));
704 if ((dummy2 = strchr (dirstr, CLOSEDIR)) != NULL)
710 if ((newd = str2dir(dirstr)) == d_invalid)
712 sprintf(errbuf, "Invalid directive %s", dirstr);
713 warning_error(wi, errbuf);
717 /* Directive found */
720 fprintf(debug, "found directive '%s'\n", dir2str(newd));
722 if (DS_Check_Order (DS, newd))
729 /* we should print here which directives should have
730 been present, and which actually are */
731 gmx_fatal(FARGS, "%s\nInvalid order for directive %s",
732 cpp_error(&handle, eCPP_SYNTAX), dir2str(newd));
738 else if (d != d_invalid)
740 /* Not a directive, just a plain string
741 * use a gigantic switch to decode,
742 * if there is a valid directive!
749 gmx_fatal(FARGS, "%s\nFound a second defaults directive.\n",
750 cpp_error(&handle, eCPP_SYNTAX));
752 bReadDefaults = TRUE;
753 nscan = sscanf(pline, "%s%s%s%lf%lf%lf",
754 nb_str, comb_str, genpairs, &fLJ, &fQQ, &fPOW);
765 get_nbparm(nb_str, comb_str, &nb_funct, &comb, wi);
766 *combination_rule = comb;
769 bGenPairs = (gmx_strncasecmp(genpairs, "Y", 1) == 0);
770 if (nb_funct != eNBF_LJ && bGenPairs)
772 gmx_fatal(FARGS, "Generating pair parameters is only supported with LJ non-bonded interactions");
788 nb_funct = ifunc_index(d_nonbond_params, nb_funct);
792 push_at(symtab, atype, batype, pline, nb_funct,
793 &nbparam, bGenPairs ? &pair : NULL, wi);
797 push_bt(d, plist, 2, NULL, batype, pline, wi);
799 case d_constrainttypes:
800 push_bt(d, plist, 2, NULL, batype, pline, wi);
805 push_nbt(d, pair, atype, pline, F_LJ14, wi);
809 push_bt(d, plist, 2, atype, NULL, pline, wi);
813 push_bt(d, plist, 3, NULL, batype, pline, wi);
815 case d_dihedraltypes:
816 /* Special routine that can read both 2 and 4 atom dihedral definitions. */
817 push_dihedraltype(d, plist, batype, pline, wi);
820 case d_nonbond_params:
821 push_nbt(d, nbparam, atype, pline, nb_funct, wi);
826 srenew(block,nblock);
827 srenew(blockinfo,nblock);
828 blk0=&(block[nblock-1]);
829 bi0=&(blockinfo[nblock-1]);
832 push_molt(symtab,bi0,pline);
836 case d_implicit_genborn_params:
837 push_gb_params(atype, pline, wi);
840 case d_implicit_surface_params:
841 gmx_fatal(FARGS, "Implicit surface directive not supported yet.");
845 push_cmaptype(d, plist, 5, atype, batype, pline, wi);
853 if (opts->couple_moltype != NULL &&
854 (opts->couple_lam0 == ecouplamNONE ||
855 opts->couple_lam0 == ecouplamQ ||
856 opts->couple_lam1 == ecouplamNONE ||
857 opts->couple_lam1 == ecouplamQ))
859 dcatt = add_atomtype_decoupled(symtab, atype,
860 &nbparam, bGenPairs ? &pair : NULL);
862 ntype = get_atomtype_ntypes(atype);
863 ncombs = (ntype*(ntype+1))/2;
864 generate_nbparams(comb, nb_funct, &(plist[nb_funct]), atype, wi);
865 ncopy = copy_nbparams(nbparam, nb_funct, &(plist[nb_funct]),
867 fprintf(stderr, "Generated %d of the %d non-bonded parameter combinations\n", ncombs-ncopy, ncombs);
868 free_nbparam(nbparam, ntype);
871 gen_pairs(&(plist[nb_funct]), &(plist[F_LJ14]), fudgeLJ, comb);
872 ncopy = copy_nbparams(pair, nb_funct, &(plist[F_LJ14]),
874 fprintf(stderr, "Generated %d of the %d 1-4 parameter combinations\n", ncombs-ncopy, ncombs);
875 free_nbparam(pair, ntype);
877 /* Copy GBSA parameters to atomtype array? */
882 push_molt(symtab, &nmol, molinfo, pline, wi);
883 srenew(block2, nmol);
884 block2[nmol-1].nr = 0;
885 mi0 = &((*molinfo)[nmol-1]);
889 push_atom(symtab, &(mi0->cgs), &(mi0->atoms), atype, pline, &lastcg, wi);
893 push_bond(d, plist, mi0->plist, &(mi0->atoms), atype, pline, FALSE,
894 bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
904 case d_position_restraints:
905 case d_angle_restraints:
906 case d_angle_restraints_z:
907 case d_distance_restraints:
908 case d_orientation_restraints:
909 case d_dihedral_restraints:
912 case d_water_polarization:
913 case d_thole_polarization:
914 push_bond(d, plist, mi0->plist, &(mi0->atoms), atype, pline, TRUE,
915 bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
918 push_cmap(d, plist, mi0->plist, &(mi0->atoms), atype, pline, wi);
922 push_vsitesn(d, mi0->plist, &(mi0->atoms), pline, wi);
926 if (!block2[nmol-1].nr)
928 init_block2(&(block2[nmol-1]), mi0->atoms.nr);
930 push_excl(pline, &(block2[nmol-1]));
934 title = put_symtab(symtab, pline);
941 push_mol(nmol, *molinfo, pline, &whichmol, &nrcopies, wi);
942 mi0 = &((*molinfo)[whichmol]);
943 srenew(molb, nmolb+1);
944 molb[nmolb].type = whichmol;
945 molb[nmolb].nmol = nrcopies;
948 bCouple = (opts->couple_moltype != NULL &&
949 (gmx_strcasecmp("system", opts->couple_moltype) == 0 ||
950 gmx_strcasecmp(*(mi0->name), opts->couple_moltype) == 0));
953 nmol_couple += nrcopies;
956 if (mi0->atoms.nr == 0)
958 gmx_fatal(FARGS, "Molecule type '%s' contains no atoms",
962 "Excluding %d bonded neighbours molecule type '%s'\n",
963 mi0->nrexcl, *mi0->name);
964 sum_q(&mi0->atoms, nrcopies, &qt, &qBt);
965 if (!mi0->bProcessed)
968 generate_excl(mi0->nrexcl,
973 merge_excl(&(mi0->excls), &(block2[whichmol]));
974 done_block2(&(block2[whichmol]));
975 make_shake(mi0->plist, &mi0->atoms, opts->nshake);
979 /* nnb contains information about first,2nd,3rd bonded neighbors.
980 * Use this to generate GB 1-2,1-3,1-4 interactions when necessary.
982 if (bGenborn == TRUE)
984 generate_gb_exclusion_interactions(mi0, atype, &nnb);
991 convert_moltype_couple(mi0, dcatt, *fudgeQQ,
992 opts->couple_lam0, opts->couple_lam1,
994 nb_funct, &(plist[nb_funct]));
996 stupid_fill_block(&mi0->mols, mi0->atoms.nr, TRUE);
997 mi0->bProcessed = TRUE;
1002 fprintf (stderr, "case: %d\n", d);
1012 status = cpp_close_file(&handle);
1013 if (status != eCPP_OK)
1015 gmx_fatal(FARGS, cpp_error(&handle, status));
1020 gmx_fio_fclose(out);
1023 if (opts->couple_moltype)
1025 if (nmol_couple == 0)
1027 gmx_fatal(FARGS, "Did not find any molecules of type '%s' for coupling",
1028 opts->couple_moltype);
1030 fprintf(stderr, "Coupling %d copies of molecule type '%s'\n",
1031 nmol_couple, opts->couple_moltype);
1034 /* this is not very clean, but fixes core dump on empty system name */
1037 title = put_symtab(symtab, "");
1039 if (fabs(qt) > 1e-4)
1041 sprintf(warn_buf, "System has non-zero total charge: %.6f\n%s\n", qt, floating_point_arithmetic_tip);
1042 warning_note(wi, warn_buf);
1044 if (fabs(qBt) > 1e-4 && !gmx_within_tol(qBt, qt, 1e-6))
1046 sprintf(warn_buf, "State B has non-zero total charge: %.6f\n%s\n", qBt, floating_point_arithmetic_tip);
1047 warning_note(wi, warn_buf);
1050 for (i = 0; i < nmol; i++)
1052 done_block2(&(block2[i]));
1056 done_bond_atomtype(&batype);
1066 char **do_top(gmx_bool bVerbose,
1067 const char *topfile,
1068 const char *topppfile,
1073 int *combination_rule,
1074 double *repulsion_power,
1076 gpp_atomtype_t atype,
1078 t_molinfo **molinfo,
1081 gmx_molblock_t **molblock,
1085 /* Tmpfile might contain a long path */
1086 const char *tmpfile;
1091 tmpfile = topppfile;
1100 printf("processing topology...\n");
1102 title = read_topol(topfile, tmpfile, opts->define, opts->include,
1103 symtab, atype, nrmols, molinfo,
1104 plist, combination_rule, repulsion_power,
1105 opts, fudgeQQ, nmolblock, molblock,
1106 ir->efep != efepNO, bGenborn, bZero, wi);
1107 if ((*combination_rule != eCOMB_GEOMETRIC) &&
1108 (ir->vdwtype == evdwUSER))
1110 warning(wi, "Using sigma/epsilon based combination rules with"
1111 " user supplied potential function may produce unwanted"
1119 static void generate_qmexcl_moltype(gmx_moltype_t *molt, unsigned char *grpnr,
1122 /* This routine expects molt->ilist to be of size F_NRE and ordered. */
1124 /* generates the exclusions between the individual QM atoms, as
1125 * these interactions should be handled by the QM subroutines and
1126 * not by the gromacs routines
1129 i, j, l, k = 0, jmax, qm_max = 0, qm_nr = 0, nratoms = 0, link_nr = 0, link_max = 0;
1131 *qm_arr = NULL, *link_arr = NULL, a1, a2, a3, a4, ftype = 0;
1137 *bQMMM, *blink, bexcl;
1139 /* First we search and select the QM atoms in an qm_arr array that
1140 * we use to create the exclusions.
1142 * we take the possibility into account that a user has defined more
1143 * than one QM group:
1145 * for that we also need to do this an ugly work-about just in case
1146 * the QM group contains the entire system...
1148 jmax = ir->opts.ngQM;
1150 /* we first search for all the QM atoms and put them in an array
1152 for (j = 0; j < jmax; j++)
1154 for (i = 0; i < molt->atoms.nr; i++)
1156 if (qm_nr >= qm_max)
1159 srenew(qm_arr, qm_max);
1161 if ((grpnr ? grpnr[i] : 0) == j)
1163 qm_arr[qm_nr++] = i;
1167 /* bQMMM[..] is an array containin TRUE/FALSE for atoms that are
1168 * QM/not QM. We first set all elements to false. Afterwards we use
1169 * the qm_arr to change the elements corresponding to the QM atoms
1172 snew(bQMMM, molt->atoms.nr);
1173 for (i = 0; i < molt->atoms.nr; i++)
1177 for (i = 0; i < qm_nr; i++)
1179 bQMMM[qm_arr[i]] = TRUE;
1182 /* We remove all bonded interactions (i.e. bonds,
1183 * angles, dihedrals, 1-4's), involving the QM atoms. The way they
1184 * are removed is as follows: if the interaction invloves 2 atoms,
1185 * it is removed if both atoms are QMatoms. If it involves 3 atoms,
1186 * it is removed if at least two of the atoms are QM atoms, if the
1187 * interaction involves 4 atoms, it is removed if there are at least
1188 * 2 QM atoms. Since this routine is called once before any forces
1189 * are computed, the top->idef.il[N].iatom[] array (see idef.h) can
1190 * be rewritten at this poitn without any problem. 25-9-2002 */
1192 /* first check weter we already have CONNBONDS: */
1193 if (molt->ilist[F_CONNBONDS].nr != 0)
1195 fprintf(stderr, "nr. of CONNBONDS present already: %d\n",
1196 molt->ilist[F_CONNBONDS].nr/3);
1197 ftype = molt->ilist[F_CONNBONDS].iatoms[0];
1198 k = molt->ilist[F_CONNBONDS].nr;
1200 /* now we delete all bonded interactions, except the ones describing
1201 * a chemical bond. These are converted to CONNBONDS
1203 for (i = 0; i < F_LJ; i++)
1205 if (i == F_CONNBONDS)
1209 nratoms = interaction_function[i].nratoms;
1211 while (j < molt->ilist[i].nr)
1217 a1 = molt->ilist[i].iatoms[j+1];
1218 a2 = molt->ilist[i].iatoms[j+2];
1219 bexcl = (bQMMM[a1] && bQMMM[a2]);
1220 /* a bonded beteen two QM atoms will be copied to the
1221 * CONNBONDS list, for reasons mentioned above
1223 if (bexcl && i < F_ANGLES)
1225 srenew(molt->ilist[F_CONNBONDS].iatoms, k+3);
1226 molt->ilist[F_CONNBONDS].nr += 3;
1227 molt->ilist[F_CONNBONDS].iatoms[k++] = ftype;
1228 molt->ilist[F_CONNBONDS].iatoms[k++] = a1;
1229 molt->ilist[F_CONNBONDS].iatoms[k++] = a2;
1233 a1 = molt->ilist[i].iatoms[j+1];
1234 a2 = molt->ilist[i].iatoms[j+2];
1235 a3 = molt->ilist[i].iatoms[j+3];
1236 bexcl = ((bQMMM[a1] && bQMMM[a2]) ||
1237 (bQMMM[a1] && bQMMM[a3]) ||
1238 (bQMMM[a2] && bQMMM[a3]));
1241 a1 = molt->ilist[i].iatoms[j+1];
1242 a2 = molt->ilist[i].iatoms[j+2];
1243 a3 = molt->ilist[i].iatoms[j+3];
1244 a4 = molt->ilist[i].iatoms[j+4];
1245 bexcl = ((bQMMM[a1] && bQMMM[a2] && bQMMM[a3]) ||
1246 (bQMMM[a1] && bQMMM[a2] && bQMMM[a4]) ||
1247 (bQMMM[a1] && bQMMM[a3] && bQMMM[a4]) ||
1248 (bQMMM[a2] && bQMMM[a3] && bQMMM[a4]));
1251 gmx_fatal(FARGS, "no such bonded interactions with %d atoms\n", nratoms);
1255 /* since the interaction involves QM atoms, these should be
1256 * removed from the MM ilist
1258 molt->ilist[i].nr -= (nratoms+1);
1259 for (l = j; l < molt->ilist[i].nr; l++)
1261 molt->ilist[i].iatoms[l] = molt->ilist[i].iatoms[l+(nratoms+1)];
1266 j += nratoms+1; /* the +1 is for the functype */
1270 /* Now, we search for atoms bonded to a QM atom because we also want
1271 * to exclude their nonbonded interactions with the QM atoms. The
1272 * reason for this is that this interaction is accounted for in the
1273 * linkatoms interaction with the QMatoms and would be counted
1276 for (i = 0; i < F_NRE; i++)
1281 while (j < molt->ilist[i].nr)
1283 a1 = molt->ilist[i].iatoms[j+1];
1284 a2 = molt->ilist[i].iatoms[j+2];
1285 if ((bQMMM[a1] && !bQMMM[a2]) || (!bQMMM[a1] && bQMMM[a2]))
1287 if (link_nr >= link_max)
1290 srenew(link_arr, link_max);
1294 link_arr[link_nr++] = a2;
1298 link_arr[link_nr++] = a1;
1305 snew(blink, molt->atoms.nr);
1306 for (i = 0; i < molt->atoms.nr; i++)
1310 for (i = 0; i < link_nr; i++)
1312 blink[link_arr[i]] = TRUE;
1314 /* creating the exclusion block for the QM atoms. Each QM atom has
1315 * as excluded elements all the other QMatoms (and itself).
1317 qmexcl.nr = molt->atoms.nr;
1318 qmexcl.nra = qm_nr*(qm_nr+link_nr)+link_nr*qm_nr;
1319 snew(qmexcl.index, qmexcl.nr+1);
1320 snew(qmexcl.a, qmexcl.nra);
1322 for (i = 0; i < qmexcl.nr; i++)
1324 qmexcl.index[i] = j;
1327 for (k = 0; k < qm_nr; k++)
1329 qmexcl.a[k+j] = qm_arr[k];
1331 for (k = 0; k < link_nr; k++)
1333 qmexcl.a[qm_nr+k+j] = link_arr[k];
1335 j += (qm_nr+link_nr);
1339 for (k = 0; k < qm_nr; k++)
1341 qmexcl.a[k+j] = qm_arr[k];
1346 qmexcl.index[qmexcl.nr] = j;
1348 /* and merging with the exclusions already present in sys.
1351 init_block2(&qmexcl2, molt->atoms.nr);
1352 b_to_b2(&qmexcl, &qmexcl2);
1353 merge_excl(&(molt->excls), &qmexcl2);
1354 done_block2(&qmexcl2);
1356 /* Finally, we also need to get rid of the pair interactions of the
1357 * classical atom bonded to the boundary QM atoms with the QMatoms,
1358 * as this interaction is already accounted for by the QM, so also
1359 * here we run the risk of double counting! We proceed in a similar
1360 * way as we did above for the other bonded interactions: */
1361 for (i = F_LJ14; i < F_COUL14; i++)
1363 nratoms = interaction_function[i].nratoms;
1365 while (j < molt->ilist[i].nr)
1368 a1 = molt->ilist[i].iatoms[j+1];
1369 a2 = molt->ilist[i].iatoms[j+2];
1370 bexcl = ((bQMMM[a1] && bQMMM[a2]) ||
1371 (blink[a1] && bQMMM[a2]) ||
1372 (bQMMM[a1] && blink[a2]));
1375 /* since the interaction involves QM atoms, these should be
1376 * removed from the MM ilist
1378 molt->ilist[i].nr -= (nratoms+1);
1379 for (k = j; k < molt->ilist[i].nr; k++)
1381 molt->ilist[i].iatoms[k] = molt->ilist[i].iatoms[k+(nratoms+1)];
1386 j += nratoms+1; /* the +1 is for the functype */
1395 } /* generate_qmexcl */
1397 void generate_qmexcl(gmx_mtop_t *sys, t_inputrec *ir, warninp_t wi)
1399 /* This routine expects molt->molt[m].ilist to be of size F_NRE and ordered.
1402 unsigned char *grpnr;
1403 int mb, mol, nat_mol, i, nr_mol_with_qm_atoms = 0;
1404 gmx_molblock_t *molb;
1407 grpnr = sys->groups.grpnr[egcQMMM];
1409 for (mb = 0; mb < sys->nmolblock; mb++)
1411 molb = &sys->molblock[mb];
1412 nat_mol = sys->moltype[molb->type].atoms.nr;
1413 for (mol = 0; mol < molb->nmol; mol++)
1416 for (i = 0; i < nat_mol; i++)
1418 if ((grpnr ? grpnr[i] : 0) < ir->opts.ngQM)
1425 nr_mol_with_qm_atoms++;
1428 /* We need to split this molblock */
1431 /* Split the molblock at this molecule */
1433 srenew(sys->molblock, sys->nmolblock);
1434 for (i = sys->nmolblock-2; i >= mb; i--)
1436 sys->molblock[i+1] = sys->molblock[i];
1438 sys->molblock[mb ].nmol = mol;
1439 sys->molblock[mb+1].nmol -= mol;
1441 molb = &sys->molblock[mb];
1445 /* Split the molblock after this molecule */
1447 srenew(sys->molblock, sys->nmolblock);
1448 molb = &sys->molblock[mb];
1449 for (i = sys->nmolblock-2; i >= mb; i--)
1451 sys->molblock[i+1] = sys->molblock[i];
1453 sys->molblock[mb ].nmol = 1;
1454 sys->molblock[mb+1].nmol -= 1;
1457 /* Add a moltype for the QMMM molecule */
1459 srenew(sys->moltype, sys->nmoltype);
1460 /* Copy the moltype struct */
1461 sys->moltype[sys->nmoltype-1] = sys->moltype[molb->type];
1462 /* Copy the exclusions to a new array, since this is the only
1463 * thing that needs to be modified for QMMM.
1465 copy_blocka(&sys->moltype[molb->type ].excls,
1466 &sys->moltype[sys->nmoltype-1].excls);
1467 /* Set the molecule type for the QMMM molblock */
1468 molb->type = sys->nmoltype - 1;
1470 generate_qmexcl_moltype(&sys->moltype[molb->type], grpnr, ir);
1478 if (nr_mol_with_qm_atoms > 1)
1480 /* generate a warning is there are QM atoms in different
1481 * topologies. In this case it is not possible at this stage to
1482 * mutualy exclude the non-bonded interactions via the
1483 * exclusions (AFAIK). Instead, the user is advised to use the
1484 * energy group exclusions in the mdp file
1487 "\nThe QM subsystem is divided over multiple topologies. "
1488 "The mutual non-bonded interactions cannot be excluded. "
1489 "There are two ways to achieve this:\n\n"
1490 "1) merge the topologies, such that the atoms of the QM "
1491 "subsystem are all present in one single topology file. "
1492 "In this case this warning will dissappear\n\n"
1493 "2) exclude the non-bonded interactions explicitly via the "
1494 "energygrp-excl option in the mdp file. if this is the case "
1495 "this warning may be ignored"