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.
49 #include <sys/types.h>
51 #include "gromacs/utility/futil.h"
53 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/fileio/gmxfio.h"
57 #include "gromacs/math/units.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/topology/block.h"
62 #include "gromacs/topology/symtab.h"
63 #include "gromacs/topology/topology.h"
64 #include "gromacs/utility/fatalerror.h"
66 #include "vsite_parm.h"
68 #include "grompp-impl.h"
72 #include "gpp_nextnb.h"
76 #include "gpp_bond_atomtype.h"
78 #include "gromacs/math/utilities.h"
80 #define OPENDIR '[' /* starting sign for directive */
81 #define CLOSEDIR ']' /* ending sign for directive */
83 static void free_nbparam(t_nbparam **param, int nr)
87 for (i = 0; i < nr; i++)
94 static int copy_nbparams(t_nbparam **param, int ftype, t_params *plist, int nr)
102 for (i = 0; i < nr; i++)
104 for (j = 0; j <= i; j++)
106 if (param[i][j].bSet)
108 for (f = 0; f < nrfp; f++)
110 plist->param[nr*i+j].c[f] = param[i][j].c[f];
111 plist->param[nr*j+i].c[f] = param[i][j].c[f];
121 static void gen_pairs(t_params *nbs, t_params *pairs, real fudge, int comb)
123 int i, j, ntp, nrfp, nrfpA, nrfpB, nnn;
128 nrfpA = interaction_function[F_LJ14].nrfpA;
129 nrfpB = interaction_function[F_LJ14].nrfpB;
132 if ((nrfp != nrfpA) || (nrfpA != nrfpB))
134 gmx_incons("Number of force parameters in gen_pairs wrong");
137 fprintf(stderr, "Generating 1-4 interactions: fudge = %g\n", fudge);
140 fprintf(debug, "Fudge factor for 1-4 interactions: %g\n", fudge);
141 fprintf(debug, "Holy Cow! there are %d types\n", ntp);
143 snew(pairs->param, pairs->nr);
144 for (i = 0; (i < ntp); i++)
147 pairs->param[i].a[0] = i / nnn;
148 pairs->param[i].a[1] = i % nnn;
149 /* Copy normal and FEP parameters and multiply by fudge factor */
153 for (j = 0; (j < nrfp); j++)
155 /* If we are using sigma/epsilon values, only the epsilon values
156 * should be scaled, but not sigma.
157 * The sigma values have even indices 0,2, etc.
159 if ((comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) && (j%2 == 0))
168 pairs->param[i].c[j] = scaling*nbs->param[i].c[j];
169 pairs->param[i].c[nrfp+j] = scaling*nbs->param[i].c[j];
174 double check_mol(gmx_mtop_t *mtop, warninp_t wi)
177 int i, mb, nmol, ri, pt;
182 /* Check mass and charge */
185 for (mb = 0; mb < mtop->nmoltype; mb++)
187 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
188 nmol = mtop->molblock[mb].nmol;
189 for (i = 0; (i < atoms->nr); i++)
191 q += nmol*atoms->atom[i].q;
192 m = atoms->atom[i].m;
193 mB = atoms->atom[i].mB;
194 pt = atoms->atom[i].ptype;
195 /* If the particle is an atom or a nucleus it must have a mass,
196 * else, if it is a shell, a vsite or a bondshell it can have mass zero
198 if (((m <= 0.0) || (mB <= 0.0)) && ((pt == eptAtom) || (pt == eptNucleus)))
200 ri = atoms->atom[i].resind;
201 sprintf(buf, "atom %s (Res %s-%d) has mass %g (state A) / %g (state B)\n",
202 *(atoms->atomname[i]),
203 *(atoms->resinfo[ri].name),
204 atoms->resinfo[ri].nr,
206 warning_error(wi, buf);
209 if (((m != 0) || (mB != 0)) && (pt == eptVSite))
211 ri = atoms->atom[i].resind;
212 sprintf(buf, "virtual site %s (Res %s-%d) has non-zero mass %g (state A) / %g (state B)\n"
213 " Check your topology.\n",
214 *(atoms->atomname[i]),
215 *(atoms->resinfo[ri].name),
216 atoms->resinfo[ri].nr,
218 warning_error(wi, buf);
219 /* The following statements make LINCS break! */
220 /* atoms->atom[i].m=0; */
227 static void sum_q(t_atoms *atoms, int n, double *qt, double *qBt)
235 for (i = 0; i < atoms->nr; i++)
237 qmolA += atoms->atom[i].q;
238 qmolB += atoms->atom[i].qB;
240 /* Unfortunately an absolute comparison,
241 * but this avoids unnecessary warnings and gmx-users mails.
243 if (fabs(qmolA) >= 1e-6 || fabs(qmolB) >= 1e-6)
250 static void get_nbparm(char *nb_str, char *comb_str, int *nb, int *comb,
254 char warn_buf[STRLEN];
257 for (i = 1; (i < eNBF_NR); i++)
259 if (gmx_strcasecmp(nb_str, enbf_names[i]) == 0)
266 *nb = strtol(nb_str, NULL, 10);
268 if ((*nb < 1) || (*nb >= eNBF_NR))
270 sprintf(warn_buf, "Invalid nonbond function selector '%s' using %s",
271 nb_str, enbf_names[1]);
272 warning_error(wi, warn_buf);
276 for (i = 1; (i < eCOMB_NR); i++)
278 if (gmx_strcasecmp(comb_str, ecomb_names[i]) == 0)
285 *comb = strtol(comb_str, NULL, 10);
287 if ((*comb < 1) || (*comb >= eCOMB_NR))
289 sprintf(warn_buf, "Invalid combination rule selector '%s' using %s",
290 comb_str, ecomb_names[1]);
291 warning_error(wi, warn_buf);
296 static char ** cpp_opts(const char *define, const char *include,
301 const char *cppadds[2];
302 char **cppopts = NULL;
303 const char *option[2] = { "-D", "-I" };
304 const char *nopt[2] = { "define", "include" };
308 char warn_buf[STRLEN];
311 cppadds[1] = include;
312 for (n = 0; (n < 2); n++)
319 while ((*ptr != '\0') && isspace(*ptr))
324 while ((*rptr != '\0') && !isspace(*rptr))
332 strncpy(buf, ptr, len);
333 if (strstr(ptr, option[n]) != ptr)
335 set_warning_line(wi, "mdp file", -1);
336 sprintf(warn_buf, "Malformed %s option %s", nopt[n], buf);
337 warning(wi, warn_buf);
341 srenew(cppopts, ++ncppopts);
342 cppopts[ncppopts-1] = strdup(buf);
350 srenew(cppopts, ++ncppopts);
351 cppopts[ncppopts-1] = NULL;
358 find_gb_bondlength(t_params *plist, int ai, int aj, real *length)
365 for (i = 0; i < F_NRE && !found; i++)
369 for (j = 0; j < plist[i].nr; j++)
371 a1 = plist[i].param[j].a[0];
372 a2 = plist[i].param[j].a[1];
374 if ( (a1 == ai && a2 == aj) || (a1 == aj && a2 == ai))
376 /* Equilibrium bond distance */
377 *length = plist[i].param[j].c[0];
390 find_gb_anglelength(t_params *plist, int ai, int ak, real *length)
392 int i, j, a1, a2, a3;
395 int status, status1, status2;
399 for (i = 0; i < F_NRE && !found; i++)
403 for (j = 0; j < plist[i].nr; j++)
405 a1 = plist[i].param[j].a[0];
406 a2 = plist[i].param[j].a[1];
407 a3 = plist[i].param[j].a[2];
409 /* We dont care what the middle atom is, but use it below */
410 if ( (a1 == ai && a3 == ak) || (a1 == ak && a3 == ai) )
412 /* Equilibrium bond distance */
413 a123 = plist[i].param[j].c[0];
414 /* Use middle atom to find reference distances r12 and r23 */
415 status1 = find_gb_bondlength(plist, a1, a2, &r12);
416 status2 = find_gb_bondlength(plist, a2, a3, &r23);
418 if (status1 == 0 && status2 == 0)
420 /* cosine theorem to get r13 */
421 *length = sqrt(r12*r12+r23*r23-(2*r12*r23*cos(a123/RAD2DEG)));
434 generate_gb_exclusion_interactions(t_molinfo *mi, gpp_atomtype_t atype, t_nextnb *nnb)
436 int i, j, k, n, ai, aj, ti, tj;
442 real radiusi, radiusj;
443 real gb_radiusi, gb_radiusj;
444 real param_c2, param_c4;
450 for (n = 1; n <= nnb->nrex; n++)
465 /* Put all higher-order exclusions into 1,4 list so we dont miss them */
472 for (ai = 0; ai < nnb->nr; ai++)
474 ti = at->atom[ai].type;
475 radiusi = get_atomtype_radius(ti, atype);
476 gb_radiusi = get_atomtype_gb_radius(ti, atype);
478 for (j = 0; j < nnb->nrexcl[ai][n]; j++)
480 aj = nnb->a[ai][n][j];
482 /* Only add the interactions once */
485 tj = at->atom[aj].type;
486 radiusj = get_atomtype_radius(tj, atype);
487 gb_radiusj = get_atomtype_gb_radius(tj, atype);
489 /* There is an exclusion of type "ftype" between atoms ai and aj */
493 /* Reference distance, not used for 1-4 interactions */
497 if (find_gb_bondlength(plist, ai, aj, &distance) != 0)
499 gmx_fatal(FARGS, "Cannot find bond length for atoms %d-%d", ai, aj);
503 if (find_gb_anglelength(plist, ai, aj, &distance) != 0)
505 gmx_fatal(FARGS, "Cannot find length for atoms %d-%d involved in angle", ai, aj);
512 /* Assign GB parameters */
514 param.c[0] = radiusi+radiusj;
515 /* Reference distance distance */
516 param.c[1] = distance;
517 /* Still parameter */
518 param.c[2] = param_c2;
520 param.c[3] = gb_radiusi+gb_radiusj;
522 param.c[4] = param_c4;
524 /* Add it to the parameter list */
525 add_param_to_list(&plist[ftype], ¶m);
536 static char **read_topol(const char *infile, const char *outfile,
537 const char *define, const char *include,
539 gpp_atomtype_t atype,
543 int *combination_rule,
548 gmx_molblock_t **molblock,
555 int i, sl, nb_funct, comb;
556 char *pline = NULL, **title = NULL;
557 char line[STRLEN], errbuf[256], comb_str[256], nb_str[256];
559 char *dirstr, *dummy2;
560 int nrcopies, nmol, nmolb = 0, nscan, ncombs, ncopy;
561 double fLJ, fQQ, fPOW;
562 gmx_molblock_t *molb = NULL;
563 t_topology *block = NULL;
564 t_molinfo *mi0 = NULL;
567 t_nbparam **nbparam, **pair;
569 real fudgeLJ = -1; /* Multiplication factor to generate 1-4 from LJ */
570 gmx_bool bReadDefaults, bReadMolType, bGenPairs, bWarn_copy_A_B;
571 double qt = 0, qBt = 0; /* total charge */
572 t_bond_atomtype batype;
574 int dcatt = -1, nmol_couple;
575 /* File handling variables */
578 char *tmp_line = NULL;
579 char warn_buf[STRLEN];
580 const char *floating_point_arithmetic_tip =
581 "Total charge should normally be an integer. See\n"
582 "http://www.gromacs.org/Documentation/Floating_Point_Arithmetic\n"
583 "for discussion on how close it should be to an integer.\n";
584 /* We need to open the output file before opening the input file,
585 * because cpp_open_file can change the current working directory.
589 out = gmx_fio_fopen(outfile, "w");
596 /* open input file */
597 status = cpp_open_file(infile, &handle, cpp_opts(define, include, wi));
600 gmx_fatal(FARGS, cpp_error(&handle, status));
603 /* some local variables */
604 DS_Init(&DS); /* directive stack */
605 nmol = 0; /* no molecules yet... */
606 d = d_invalid; /* first thing should be a directive */
607 nbparam = NULL; /* The temporary non-bonded matrix */
608 pair = NULL; /* The temporary pair interaction matrix */
609 block2 = NULL; /* the extra exclusions */
611 *reppow = 12.0; /* Default value for repulsion power */
615 /* Init the number of CMAP torsion angles and grid spacing */
616 plist->grid_spacing = 0;
619 bWarn_copy_A_B = bFEP;
621 batype = init_bond_atomtype();
622 /* parse the actual file */
623 bReadDefaults = FALSE;
625 bReadMolType = FALSE;
630 status = cpp_read_line(&handle, STRLEN, line);
631 done = (status == eCPP_EOF);
634 if (status != eCPP_OK)
636 gmx_fatal(FARGS, cpp_error(&handle, status));
640 fprintf(out, "%s\n", line);
643 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
645 pline = strdup(line);
647 /* Strip trailing '\' from pline, if it exists */
649 if ((sl > 0) && (pline[sl-1] == CONTINUE))
654 /* build one long line from several fragments - necessary for CMAP */
655 while (continuing(line))
657 status = cpp_read_line(&handle, STRLEN, line);
658 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
660 /* Since we depend on the '\' being present to continue to read, we copy line
661 * to a tmp string, strip the '\' from that string, and cat it to pline
663 tmp_line = strdup(line);
665 sl = strlen(tmp_line);
666 if ((sl > 0) && (tmp_line[sl-1] == CONTINUE))
668 tmp_line[sl-1] = ' ';
671 done = (status == eCPP_EOF);
674 if (status != eCPP_OK)
676 gmx_fatal(FARGS, cpp_error(&handle, status));
680 fprintf(out, "%s\n", line);
684 srenew(pline, strlen(pline)+strlen(tmp_line)+1);
685 strcat(pline, tmp_line);
689 /* skip trailing and leading spaces and comment text */
690 strip_comment (pline);
693 /* if there is something left... */
694 if ((int)strlen(pline) > 0)
696 if (pline[0] == OPENDIR)
698 /* A directive on this line: copy the directive
699 * without the brackets into dirstr, then
700 * skip spaces and tabs on either side of directive
702 dirstr = strdup((pline+1));
703 if ((dummy2 = strchr (dirstr, CLOSEDIR)) != NULL)
709 if ((newd = str2dir(dirstr)) == d_invalid)
711 sprintf(errbuf, "Invalid directive %s", dirstr);
712 warning_error(wi, errbuf);
716 /* Directive found */
719 fprintf(debug, "found directive '%s'\n", dir2str(newd));
721 if (DS_Check_Order (DS, newd))
728 /* we should print here which directives should have
729 been present, and which actually are */
730 gmx_fatal(FARGS, "%s\nInvalid order for directive %s",
731 cpp_error(&handle, eCPP_SYNTAX), dir2str(newd));
737 else if (d != d_invalid)
739 /* Not a directive, just a plain string
740 * use a gigantic switch to decode,
741 * if there is a valid directive!
748 gmx_fatal(FARGS, "%s\nFound a second defaults directive.\n",
749 cpp_error(&handle, eCPP_SYNTAX));
751 bReadDefaults = TRUE;
752 nscan = sscanf(pline, "%s%s%s%lf%lf%lf",
753 nb_str, comb_str, genpairs, &fLJ, &fQQ, &fPOW);
764 get_nbparm(nb_str, comb_str, &nb_funct, &comb, wi);
765 *combination_rule = comb;
768 bGenPairs = (gmx_strncasecmp(genpairs, "Y", 1) == 0);
769 if (nb_funct != eNBF_LJ && bGenPairs)
771 gmx_fatal(FARGS, "Generating pair parameters is only supported with LJ non-bonded interactions");
787 nb_funct = ifunc_index(d_nonbond_params, nb_funct);
791 push_at(symtab, atype, batype, pline, nb_funct,
792 &nbparam, bGenPairs ? &pair : NULL, wi);
796 push_bt(d, plist, 2, NULL, batype, pline, wi);
798 case d_constrainttypes:
799 push_bt(d, plist, 2, NULL, batype, pline, wi);
804 push_nbt(d, pair, atype, pline, F_LJ14, wi);
808 push_bt(d, plist, 2, atype, NULL, pline, wi);
812 push_bt(d, plist, 3, NULL, batype, pline, wi);
814 case d_dihedraltypes:
815 /* Special routine that can read both 2 and 4 atom dihedral definitions. */
816 push_dihedraltype(d, plist, batype, pline, wi);
819 case d_nonbond_params:
820 push_nbt(d, nbparam, atype, pline, nb_funct, wi);
825 srenew(block,nblock);
826 srenew(blockinfo,nblock);
827 blk0=&(block[nblock-1]);
828 bi0=&(blockinfo[nblock-1]);
831 push_molt(symtab,bi0,pline);
835 case d_implicit_genborn_params:
836 push_gb_params(atype, pline, wi);
839 case d_implicit_surface_params:
840 gmx_fatal(FARGS, "Implicit surface directive not supported yet.");
844 push_cmaptype(d, plist, 5, atype, batype, pline, wi);
852 if (opts->couple_moltype != NULL &&
853 (opts->couple_lam0 == ecouplamNONE ||
854 opts->couple_lam0 == ecouplamQ ||
855 opts->couple_lam1 == ecouplamNONE ||
856 opts->couple_lam1 == ecouplamQ))
858 dcatt = add_atomtype_decoupled(symtab, atype,
859 &nbparam, bGenPairs ? &pair : NULL);
861 ntype = get_atomtype_ntypes(atype);
862 ncombs = (ntype*(ntype+1))/2;
863 generate_nbparams(comb, nb_funct, &(plist[nb_funct]), atype, wi);
864 ncopy = copy_nbparams(nbparam, nb_funct, &(plist[nb_funct]),
866 fprintf(stderr, "Generated %d of the %d non-bonded parameter combinations\n", ncombs-ncopy, ncombs);
867 free_nbparam(nbparam, ntype);
870 gen_pairs(&(plist[nb_funct]), &(plist[F_LJ14]), fudgeLJ, comb);
871 ncopy = copy_nbparams(pair, nb_funct, &(plist[F_LJ14]),
873 fprintf(stderr, "Generated %d of the %d 1-4 parameter combinations\n", ncombs-ncopy, ncombs);
874 free_nbparam(pair, ntype);
876 /* Copy GBSA parameters to atomtype array? */
881 push_molt(symtab, &nmol, molinfo, pline, wi);
882 srenew(block2, nmol);
883 block2[nmol-1].nr = 0;
884 mi0 = &((*molinfo)[nmol-1]);
888 push_atom(symtab, &(mi0->cgs), &(mi0->atoms), atype, pline, &lastcg, wi);
892 push_bond(d, plist, mi0->plist, &(mi0->atoms), atype, pline, FALSE,
893 bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
903 case d_position_restraints:
904 case d_angle_restraints:
905 case d_angle_restraints_z:
906 case d_distance_restraints:
907 case d_orientation_restraints:
908 case d_dihedral_restraints:
911 case d_water_polarization:
912 case d_thole_polarization:
913 push_bond(d, plist, mi0->plist, &(mi0->atoms), atype, pline, TRUE,
914 bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
917 push_cmap(d, plist, mi0->plist, &(mi0->atoms), atype, pline, wi);
921 push_vsitesn(d, mi0->plist, &(mi0->atoms), pline, wi);
925 if (!block2[nmol-1].nr)
927 init_block2(&(block2[nmol-1]), mi0->atoms.nr);
929 push_excl(pline, &(block2[nmol-1]));
933 title = put_symtab(symtab, pline);
940 push_mol(nmol, *molinfo, pline, &whichmol, &nrcopies, wi);
941 mi0 = &((*molinfo)[whichmol]);
942 srenew(molb, nmolb+1);
943 molb[nmolb].type = whichmol;
944 molb[nmolb].nmol = nrcopies;
947 bCouple = (opts->couple_moltype != NULL &&
948 (gmx_strcasecmp("system", opts->couple_moltype) == 0 ||
949 gmx_strcasecmp(*(mi0->name), opts->couple_moltype) == 0));
952 nmol_couple += nrcopies;
955 if (mi0->atoms.nr == 0)
957 gmx_fatal(FARGS, "Molecule type '%s' contains no atoms",
961 "Excluding %d bonded neighbours molecule type '%s'\n",
962 mi0->nrexcl, *mi0->name);
963 sum_q(&mi0->atoms, nrcopies, &qt, &qBt);
964 if (!mi0->bProcessed)
967 generate_excl(mi0->nrexcl,
972 merge_excl(&(mi0->excls), &(block2[whichmol]));
973 done_block2(&(block2[whichmol]));
974 make_shake(mi0->plist, &mi0->atoms, opts->nshake);
978 /* nnb contains information about first,2nd,3rd bonded neighbors.
979 * Use this to generate GB 1-2,1-3,1-4 interactions when necessary.
981 if (bGenborn == TRUE)
983 generate_gb_exclusion_interactions(mi0, atype, &nnb);
990 convert_moltype_couple(mi0, dcatt, *fudgeQQ,
991 opts->couple_lam0, opts->couple_lam1,
993 nb_funct, &(plist[nb_funct]));
995 stupid_fill_block(&mi0->mols, mi0->atoms.nr, TRUE);
996 mi0->bProcessed = TRUE;
1001 fprintf (stderr, "case: %d\n", d);
1002 gmx_incons("unknown directive");
1011 status = cpp_close_file(&handle);
1012 if (status != eCPP_OK)
1014 gmx_fatal(FARGS, cpp_error(&handle, status));
1019 gmx_fio_fclose(out);
1022 if (opts->couple_moltype)
1024 if (nmol_couple == 0)
1026 gmx_fatal(FARGS, "Did not find any molecules of type '%s' for coupling",
1027 opts->couple_moltype);
1029 fprintf(stderr, "Coupling %d copies of molecule type '%s'\n",
1030 nmol_couple, opts->couple_moltype);
1033 /* this is not very clean, but fixes core dump on empty system name */
1036 title = put_symtab(symtab, "");
1038 if (fabs(qt) > 1e-4)
1040 sprintf(warn_buf, "System has non-zero total charge: %.6f\n%s\n", qt, floating_point_arithmetic_tip);
1041 warning_note(wi, warn_buf);
1043 if (fabs(qBt) > 1e-4 && !gmx_within_tol(qBt, qt, 1e-6))
1045 sprintf(warn_buf, "State B has non-zero total charge: %.6f\n%s\n", qBt, floating_point_arithmetic_tip);
1046 warning_note(wi, warn_buf);
1049 for (i = 0; i < nmol; i++)
1051 done_block2(&(block2[i]));
1055 done_bond_atomtype(&batype);
1065 char **do_top(gmx_bool bVerbose,
1066 const char *topfile,
1067 const char *topppfile,
1072 int *combination_rule,
1073 double *repulsion_power,
1075 gpp_atomtype_t atype,
1077 t_molinfo **molinfo,
1080 gmx_molblock_t **molblock,
1084 /* Tmpfile might contain a long path */
1085 const char *tmpfile;
1090 tmpfile = topppfile;
1099 printf("processing topology...\n");
1101 title = read_topol(topfile, tmpfile, opts->define, opts->include,
1102 symtab, atype, nrmols, molinfo,
1103 plist, combination_rule, repulsion_power,
1104 opts, fudgeQQ, nmolblock, molblock,
1105 ir->efep != efepNO, bGenborn, bZero, wi);
1106 if ((*combination_rule != eCOMB_GEOMETRIC) &&
1107 (ir->vdwtype == evdwUSER))
1109 warning(wi, "Using sigma/epsilon based combination rules with"
1110 " user supplied potential function may produce unwanted"
1118 static void generate_qmexcl_moltype(gmx_moltype_t *molt, unsigned char *grpnr,
1121 /* This routine expects molt->ilist to be of size F_NRE and ordered. */
1123 /* generates the exclusions between the individual QM atoms, as
1124 * these interactions should be handled by the QM subroutines and
1125 * not by the gromacs routines
1128 i, j, l, k = 0, jmax, qm_max = 0, qm_nr = 0, nratoms = 0, link_nr = 0, link_max = 0;
1130 *qm_arr = NULL, *link_arr = NULL, a1, a2, a3, a4, ftype = 0;
1136 *bQMMM, *blink, bexcl;
1138 /* First we search and select the QM atoms in an qm_arr array that
1139 * we use to create the exclusions.
1141 * we take the possibility into account that a user has defined more
1142 * than one QM group:
1144 * for that we also need to do this an ugly work-about just in case
1145 * the QM group contains the entire system...
1147 jmax = ir->opts.ngQM;
1149 /* we first search for all the QM atoms and put them in an array
1151 for (j = 0; j < jmax; j++)
1153 for (i = 0; i < molt->atoms.nr; i++)
1155 if (qm_nr >= qm_max)
1158 srenew(qm_arr, qm_max);
1160 if ((grpnr ? grpnr[i] : 0) == j)
1162 qm_arr[qm_nr++] = i;
1166 /* bQMMM[..] is an array containin TRUE/FALSE for atoms that are
1167 * QM/not QM. We first set all elements to false. Afterwards we use
1168 * the qm_arr to change the elements corresponding to the QM atoms
1171 snew(bQMMM, molt->atoms.nr);
1172 for (i = 0; i < molt->atoms.nr; i++)
1176 for (i = 0; i < qm_nr; i++)
1178 bQMMM[qm_arr[i]] = TRUE;
1181 /* We remove all bonded interactions (i.e. bonds,
1182 * angles, dihedrals, 1-4's), involving the QM atoms. The way they
1183 * are removed is as follows: if the interaction invloves 2 atoms,
1184 * it is removed if both atoms are QMatoms. If it involves 3 atoms,
1185 * it is removed if at least two of the atoms are QM atoms, if the
1186 * interaction involves 4 atoms, it is removed if there are at least
1187 * 2 QM atoms. Since this routine is called once before any forces
1188 * are computed, the top->idef.il[N].iatom[] array (see idef.h) can
1189 * be rewritten at this poitn without any problem. 25-9-2002 */
1191 /* first check weter we already have CONNBONDS: */
1192 if (molt->ilist[F_CONNBONDS].nr != 0)
1194 fprintf(stderr, "nr. of CONNBONDS present already: %d\n",
1195 molt->ilist[F_CONNBONDS].nr/3);
1196 ftype = molt->ilist[F_CONNBONDS].iatoms[0];
1197 k = molt->ilist[F_CONNBONDS].nr;
1199 /* now we delete all bonded interactions, except the ones describing
1200 * a chemical bond. These are converted to CONNBONDS
1202 for (i = 0; i < F_LJ; i++)
1204 if (i == F_CONNBONDS)
1208 nratoms = interaction_function[i].nratoms;
1210 while (j < molt->ilist[i].nr)
1216 a1 = molt->ilist[i].iatoms[j+1];
1217 a2 = molt->ilist[i].iatoms[j+2];
1218 bexcl = (bQMMM[a1] && bQMMM[a2]);
1219 /* a bonded beteen two QM atoms will be copied to the
1220 * CONNBONDS list, for reasons mentioned above
1222 if (bexcl && i < F_ANGLES)
1224 srenew(molt->ilist[F_CONNBONDS].iatoms, k+3);
1225 molt->ilist[F_CONNBONDS].nr += 3;
1226 molt->ilist[F_CONNBONDS].iatoms[k++] = ftype;
1227 molt->ilist[F_CONNBONDS].iatoms[k++] = a1;
1228 molt->ilist[F_CONNBONDS].iatoms[k++] = a2;
1232 a1 = molt->ilist[i].iatoms[j+1];
1233 a2 = molt->ilist[i].iatoms[j+2];
1234 a3 = molt->ilist[i].iatoms[j+3];
1235 bexcl = ((bQMMM[a1] && bQMMM[a2]) ||
1236 (bQMMM[a1] && bQMMM[a3]) ||
1237 (bQMMM[a2] && bQMMM[a3]));
1240 a1 = molt->ilist[i].iatoms[j+1];
1241 a2 = molt->ilist[i].iatoms[j+2];
1242 a3 = molt->ilist[i].iatoms[j+3];
1243 a4 = molt->ilist[i].iatoms[j+4];
1244 bexcl = ((bQMMM[a1] && bQMMM[a2] && bQMMM[a3]) ||
1245 (bQMMM[a1] && bQMMM[a2] && bQMMM[a4]) ||
1246 (bQMMM[a1] && bQMMM[a3] && bQMMM[a4]) ||
1247 (bQMMM[a2] && bQMMM[a3] && bQMMM[a4]));
1250 gmx_fatal(FARGS, "no such bonded interactions with %d atoms\n", nratoms);
1254 /* since the interaction involves QM atoms, these should be
1255 * removed from the MM ilist
1257 molt->ilist[i].nr -= (nratoms+1);
1258 for (l = j; l < molt->ilist[i].nr; l++)
1260 molt->ilist[i].iatoms[l] = molt->ilist[i].iatoms[l+(nratoms+1)];
1265 j += nratoms+1; /* the +1 is for the functype */
1269 /* Now, we search for atoms bonded to a QM atom because we also want
1270 * to exclude their nonbonded interactions with the QM atoms. The
1271 * reason for this is that this interaction is accounted for in the
1272 * linkatoms interaction with the QMatoms and would be counted
1275 for (i = 0; i < F_NRE; i++)
1280 while (j < molt->ilist[i].nr)
1282 a1 = molt->ilist[i].iatoms[j+1];
1283 a2 = molt->ilist[i].iatoms[j+2];
1284 if ((bQMMM[a1] && !bQMMM[a2]) || (!bQMMM[a1] && bQMMM[a2]))
1286 if (link_nr >= link_max)
1289 srenew(link_arr, link_max);
1293 link_arr[link_nr++] = a2;
1297 link_arr[link_nr++] = a1;
1304 snew(blink, molt->atoms.nr);
1305 for (i = 0; i < molt->atoms.nr; i++)
1309 for (i = 0; i < link_nr; i++)
1311 blink[link_arr[i]] = TRUE;
1313 /* creating the exclusion block for the QM atoms. Each QM atom has
1314 * as excluded elements all the other QMatoms (and itself).
1316 qmexcl.nr = molt->atoms.nr;
1317 qmexcl.nra = qm_nr*(qm_nr+link_nr)+link_nr*qm_nr;
1318 snew(qmexcl.index, qmexcl.nr+1);
1319 snew(qmexcl.a, qmexcl.nra);
1321 for (i = 0; i < qmexcl.nr; i++)
1323 qmexcl.index[i] = j;
1326 for (k = 0; k < qm_nr; k++)
1328 qmexcl.a[k+j] = qm_arr[k];
1330 for (k = 0; k < link_nr; k++)
1332 qmexcl.a[qm_nr+k+j] = link_arr[k];
1334 j += (qm_nr+link_nr);
1338 for (k = 0; k < qm_nr; k++)
1340 qmexcl.a[k+j] = qm_arr[k];
1345 qmexcl.index[qmexcl.nr] = j;
1347 /* and merging with the exclusions already present in sys.
1350 init_block2(&qmexcl2, molt->atoms.nr);
1351 b_to_b2(&qmexcl, &qmexcl2);
1352 merge_excl(&(molt->excls), &qmexcl2);
1353 done_block2(&qmexcl2);
1355 /* Finally, we also need to get rid of the pair interactions of the
1356 * classical atom bonded to the boundary QM atoms with the QMatoms,
1357 * as this interaction is already accounted for by the QM, so also
1358 * here we run the risk of double counting! We proceed in a similar
1359 * way as we did above for the other bonded interactions: */
1360 for (i = F_LJ14; i < F_COUL14; i++)
1362 nratoms = interaction_function[i].nratoms;
1364 while (j < molt->ilist[i].nr)
1366 a1 = molt->ilist[i].iatoms[j+1];
1367 a2 = molt->ilist[i].iatoms[j+2];
1368 bexcl = ((bQMMM[a1] && bQMMM[a2]) ||
1369 (blink[a1] && bQMMM[a2]) ||
1370 (bQMMM[a1] && blink[a2]));
1373 /* since the interaction involves QM atoms, these should be
1374 * removed from the MM ilist
1376 molt->ilist[i].nr -= (nratoms+1);
1377 for (k = j; k < molt->ilist[i].nr; k++)
1379 molt->ilist[i].iatoms[k] = molt->ilist[i].iatoms[k+(nratoms+1)];
1384 j += nratoms+1; /* the +1 is for the functype */
1393 } /* generate_qmexcl */
1395 void generate_qmexcl(gmx_mtop_t *sys, t_inputrec *ir, warninp_t wi)
1397 /* This routine expects molt->molt[m].ilist to be of size F_NRE and ordered.
1400 unsigned char *grpnr;
1401 int mb, mol, nat_mol, i, nr_mol_with_qm_atoms = 0;
1402 gmx_molblock_t *molb;
1405 grpnr = sys->groups.grpnr[egcQMMM];
1407 for (mb = 0; mb < sys->nmolblock; mb++)
1409 molb = &sys->molblock[mb];
1410 nat_mol = sys->moltype[molb->type].atoms.nr;
1411 for (mol = 0; mol < molb->nmol; mol++)
1414 for (i = 0; i < nat_mol; i++)
1416 if ((grpnr ? grpnr[i] : 0) < ir->opts.ngQM)
1423 nr_mol_with_qm_atoms++;
1426 /* We need to split this molblock */
1429 /* Split the molblock at this molecule */
1431 srenew(sys->molblock, sys->nmolblock);
1432 for (i = sys->nmolblock-2; i >= mb; i--)
1434 sys->molblock[i+1] = sys->molblock[i];
1436 sys->molblock[mb ].nmol = mol;
1437 sys->molblock[mb+1].nmol -= mol;
1439 molb = &sys->molblock[mb];
1443 /* Split the molblock after this molecule */
1445 srenew(sys->molblock, sys->nmolblock);
1446 molb = &sys->molblock[mb];
1447 for (i = sys->nmolblock-2; i >= mb; i--)
1449 sys->molblock[i+1] = sys->molblock[i];
1451 sys->molblock[mb ].nmol = 1;
1452 sys->molblock[mb+1].nmol -= 1;
1455 /* Add a moltype for the QMMM molecule */
1457 srenew(sys->moltype, sys->nmoltype);
1458 /* Copy the moltype struct */
1459 sys->moltype[sys->nmoltype-1] = sys->moltype[molb->type];
1460 /* Copy the exclusions to a new array, since this is the only
1461 * thing that needs to be modified for QMMM.
1463 copy_blocka(&sys->moltype[molb->type ].excls,
1464 &sys->moltype[sys->nmoltype-1].excls);
1465 /* Set the molecule type for the QMMM molblock */
1466 molb->type = sys->nmoltype - 1;
1468 generate_qmexcl_moltype(&sys->moltype[molb->type], grpnr, ir);
1476 if (nr_mol_with_qm_atoms > 1)
1478 /* generate a warning is there are QM atoms in different
1479 * topologies. In this case it is not possible at this stage to
1480 * mutualy exclude the non-bonded interactions via the
1481 * exclusions (AFAIK). Instead, the user is advised to use the
1482 * energy group exclusions in the mdp file
1485 "\nThe QM subsystem is divided over multiple topologies. "
1486 "The mutual non-bonded interactions cannot be excluded. "
1487 "There are two ways to achieve this:\n\n"
1488 "1) merge the topologies, such that the atoms of the QM "
1489 "subsystem are all present in one single topology file. "
1490 "In this case this warning will dissappear\n\n"
1491 "2) exclude the non-bonded interactions explicitly via the "
1492 "energygrp-excl option in the mdp file. if this is the case "
1493 "this warning may be ignored"