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 OPENDIR '[' /* starting sign for directive */
78 #define CLOSEDIR ']' /* ending sign for directive */
80 static void free_nbparam(t_nbparam **param, int nr)
84 for (i = 0; i < nr; i++)
91 static int copy_nbparams(t_nbparam **param, int ftype, t_params *plist, int nr)
99 for (i = 0; i < nr; i++)
101 for (j = 0; j <= i; j++)
103 if (param[i][j].bSet)
105 for (f = 0; f < nrfp; f++)
107 plist->param[nr*i+j].c[f] = param[i][j].c[f];
108 plist->param[nr*j+i].c[f] = param[i][j].c[f];
118 static void gen_pairs(t_params *nbs, t_params *pairs, real fudge, int comb)
120 int i, j, ntp, nrfp, nrfpA, nrfpB, nnn;
125 nrfpA = interaction_function[F_LJ14].nrfpA;
126 nrfpB = interaction_function[F_LJ14].nrfpB;
129 if ((nrfp != nrfpA) || (nrfpA != nrfpB))
131 gmx_incons("Number of force parameters in gen_pairs wrong");
134 fprintf(stderr, "Generating 1-4 interactions: fudge = %g\n", fudge);
137 fprintf(debug, "Fudge factor for 1-4 interactions: %g\n", fudge);
138 fprintf(debug, "Holy Cow! there are %d types\n", ntp);
140 snew(pairs->param, pairs->nr);
141 for (i = 0; (i < ntp); i++)
144 pairs->param[i].a[0] = i / nnn;
145 pairs->param[i].a[1] = i % nnn;
146 /* Copy normal and FEP parameters and multiply by fudge factor */
150 for (j = 0; (j < nrfp); j++)
152 /* If we are using sigma/epsilon values, only the epsilon values
153 * should be scaled, but not sigma.
154 * The sigma values have even indices 0,2, etc.
156 if ((comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) && (j%2 == 0))
165 pairs->param[i].c[j] = scaling*nbs->param[i].c[j];
166 pairs->param[i].c[nrfp+j] = scaling*nbs->param[i].c[j];
171 double check_mol(gmx_mtop_t *mtop, warninp_t wi)
174 int i, mb, nmol, ri, pt;
179 /* Check mass and charge */
182 for (mb = 0; mb < mtop->nmoltype; mb++)
184 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
185 nmol = mtop->molblock[mb].nmol;
186 for (i = 0; (i < atoms->nr); i++)
188 q += nmol*atoms->atom[i].q;
189 m = atoms->atom[i].m;
190 pt = atoms->atom[i].ptype;
191 /* If the particle is an atom or a nucleus it must have a mass,
192 * else, if it is a shell, a vsite or a bondshell it can have mass zero
194 if ((m <= 0.0) && ((pt == eptAtom) || (pt == eptNucleus)))
196 ri = atoms->atom[i].resind;
197 sprintf(buf, "atom %s (Res %s-%d) has mass %g\n",
198 *(atoms->atomname[i]),
199 *(atoms->resinfo[ri].name),
200 atoms->resinfo[ri].nr,
202 warning_error(wi, buf);
205 if ((m != 0) && (pt == eptVSite))
207 ri = atoms->atom[i].resind;
208 sprintf(buf, "virtual site %s (Res %s-%d) has non-zero mass %g\n"
209 " Check your topology.\n",
210 *(atoms->atomname[i]),
211 *(atoms->resinfo[ri].name),
212 atoms->resinfo[ri].nr,
214 warning_error(wi, buf);
215 /* The following statements make LINCS break! */
216 /* atoms->atom[i].m=0; */
223 static void sum_q(t_atoms *atoms, int n, double *qt, double *qBt)
231 for (i = 0; i < atoms->nr; i++)
233 qmolA += atoms->atom[i].q;
234 qmolB += atoms->atom[i].qB;
236 /* Unfortunately an absolute comparison,
237 * but this avoids unnecessary warnings and gmx-users mails.
239 if (fabs(qmolA) >= 1e-6 || fabs(qmolB) >= 1e-6)
246 static void get_nbparm(char *nb_str, char *comb_str, int *nb, int *comb,
250 char warn_buf[STRLEN];
253 for (i = 1; (i < eNBF_NR); i++)
255 if (gmx_strcasecmp(nb_str, enbf_names[i]) == 0)
262 *nb = strtol(nb_str, NULL, 10);
264 if ((*nb < 1) || (*nb >= eNBF_NR))
266 sprintf(warn_buf, "Invalid nonbond function selector '%s' using %s",
267 nb_str, enbf_names[1]);
268 warning_error(wi, warn_buf);
272 for (i = 1; (i < eCOMB_NR); i++)
274 if (gmx_strcasecmp(comb_str, ecomb_names[i]) == 0)
281 *comb = strtol(comb_str, NULL, 10);
283 if ((*comb < 1) || (*comb >= eCOMB_NR))
285 sprintf(warn_buf, "Invalid combination rule selector '%s' using %s",
286 comb_str, ecomb_names[1]);
287 warning_error(wi, warn_buf);
292 static char ** cpp_opts(const char *define, const char *include,
297 const char *cppadds[2];
298 char **cppopts = NULL;
299 const char *option[2] = { "-D", "-I" };
300 const char *nopt[2] = { "define", "include" };
304 char warn_buf[STRLEN];
307 cppadds[1] = include;
308 for (n = 0; (n < 2); n++)
315 while ((*ptr != '\0') && isspace(*ptr))
320 while ((*rptr != '\0') && !isspace(*rptr))
328 strncpy(buf, ptr, len);
329 if (strstr(ptr, option[n]) != ptr)
331 set_warning_line(wi, "mdp file", -1);
332 sprintf(warn_buf, "Malformed %s option %s", nopt[n], buf);
333 warning(wi, warn_buf);
337 srenew(cppopts, ++ncppopts);
338 cppopts[ncppopts-1] = strdup(buf);
346 srenew(cppopts, ++ncppopts);
347 cppopts[ncppopts-1] = NULL;
354 find_gb_bondlength(t_params *plist, int ai, int aj, real *length)
361 for (i = 0; i < F_NRE && !found; i++)
365 for (j = 0; j < plist[i].nr; j++)
367 a1 = plist[i].param[j].a[0];
368 a2 = plist[i].param[j].a[1];
370 if ( (a1 == ai && a2 == aj) || (a1 == aj && a2 == ai))
372 /* Equilibrium bond distance */
373 *length = plist[i].param[j].c[0];
386 find_gb_anglelength(t_params *plist, int ai, int ak, real *length)
388 int i, j, a1, a2, a3;
391 int status, status1, status2;
395 for (i = 0; i < F_NRE && !found; i++)
399 for (j = 0; j < plist[i].nr; j++)
401 a1 = plist[i].param[j].a[0];
402 a2 = plist[i].param[j].a[1];
403 a3 = plist[i].param[j].a[2];
405 /* We dont care what the middle atom is, but use it below */
406 if ( (a1 == ai && a3 == ak) || (a1 == ak && a3 == ai) )
408 /* Equilibrium bond distance */
409 a123 = plist[i].param[j].c[0];
410 /* Use middle atom to find reference distances r12 and r23 */
411 status1 = find_gb_bondlength(plist, a1, a2, &r12);
412 status2 = find_gb_bondlength(plist, a2, a3, &r23);
414 if (status1 == 0 && status2 == 0)
416 /* cosine theorem to get r13 */
417 *length = sqrt(r12*r12+r23*r23-(2*r12*r23*cos(a123/RAD2DEG)));
430 generate_gb_exclusion_interactions(t_molinfo *mi, gpp_atomtype_t atype, t_nextnb *nnb)
432 int i, j, k, n, ai, aj, ti, tj;
438 real radiusi, radiusj;
439 real gb_radiusi, gb_radiusj;
440 real param_c2, param_c4;
446 for (n = 1; n <= nnb->nrex; n++)
461 /* Put all higher-order exclusions into 1,4 list so we dont miss them */
468 for (ai = 0; ai < nnb->nr; ai++)
470 ti = at->atom[ai].type;
471 radiusi = get_atomtype_radius(ti, atype);
472 gb_radiusi = get_atomtype_gb_radius(ti, atype);
474 for (j = 0; j < nnb->nrexcl[ai][n]; j++)
476 aj = nnb->a[ai][n][j];
478 /* Only add the interactions once */
481 tj = at->atom[aj].type;
482 radiusj = get_atomtype_radius(tj, atype);
483 gb_radiusj = get_atomtype_gb_radius(tj, atype);
485 /* There is an exclusion of type "ftype" between atoms ai and aj */
489 /* Reference distance, not used for 1-4 interactions */
493 if (find_gb_bondlength(plist, ai, aj, &distance) != 0)
495 gmx_fatal(FARGS, "Cannot find bond length for atoms %d-%d", ai, aj);
499 if (find_gb_anglelength(plist, ai, aj, &distance) != 0)
501 gmx_fatal(FARGS, "Cannot find length for atoms %d-%d involved in angle", ai, aj);
508 /* Assign GB parameters */
510 param.c[0] = radiusi+radiusj;
511 /* Reference distance distance */
512 param.c[1] = distance;
513 /* Still parameter */
514 param.c[2] = param_c2;
516 param.c[3] = gb_radiusi+gb_radiusj;
518 param.c[4] = param_c4;
520 /* Add it to the parameter list */
521 add_param_to_list(&plist[ftype], ¶m);
532 static char **read_topol(const char *infile, const char *outfile,
533 const char *define, const char *include,
535 gpp_atomtype_t atype,
539 int *combination_rule,
544 gmx_molblock_t **molblock,
551 int i, sl, nb_funct, comb;
552 char *pline = NULL, **title = NULL;
553 char line[STRLEN], errbuf[256], comb_str[256], nb_str[256];
555 char *dirstr, *dummy2;
556 int nrcopies, nmol, nmolb = 0, nscan, ncombs, ncopy;
557 double fLJ, fQQ, fPOW;
558 gmx_molblock_t *molb = NULL;
559 t_topology *block = NULL;
560 t_molinfo *mi0 = NULL;
563 t_nbparam **nbparam, **pair;
565 real fudgeLJ = -1; /* Multiplication factor to generate 1-4 from LJ */
566 gmx_bool bReadDefaults, bReadMolType, bGenPairs, bWarn_copy_A_B;
567 double qt = 0, qBt = 0; /* total charge */
568 t_bond_atomtype batype;
570 int dcatt = -1, nmol_couple;
571 /* File handling variables */
574 char *tmp_line = NULL;
575 char warn_buf[STRLEN];
576 const char *floating_point_arithmetic_tip =
577 "Total charge should normally be an integer. See\n"
578 "http://www.gromacs.org/Documentation/Floating_Point_Arithmetic\n"
579 "for discussion on how close it should be to an integer.\n";
580 /* We need to open the output file before opening the input file,
581 * because cpp_open_file can change the current working directory.
585 out = gmx_fio_fopen(outfile, "w");
592 /* open input file */
593 status = cpp_open_file(infile, &handle, cpp_opts(define, include, wi));
596 gmx_fatal(FARGS, cpp_error(&handle, status));
599 /* some local variables */
600 DS_Init(&DS); /* directive stack */
601 nmol = 0; /* no molecules yet... */
602 d = d_invalid; /* first thing should be a directive */
603 nbparam = NULL; /* The temporary non-bonded matrix */
604 pair = NULL; /* The temporary pair interaction matrix */
605 block2 = NULL; /* the extra exclusions */
607 *reppow = 12.0; /* Default value for repulsion power */
611 /* Init the number of CMAP torsion angles and grid spacing */
612 plist->grid_spacing = 0;
615 bWarn_copy_A_B = bFEP;
617 batype = init_bond_atomtype();
618 /* parse the actual file */
619 bReadDefaults = FALSE;
621 bReadMolType = FALSE;
626 status = cpp_read_line(&handle, STRLEN, line);
627 done = (status == eCPP_EOF);
630 if (status != eCPP_OK)
632 gmx_fatal(FARGS, cpp_error(&handle, status));
636 fprintf(out, "%s\n", line);
639 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
641 pline = strdup(line);
643 /* Strip trailing '\' from pline, if it exists */
645 if ((sl > 0) && (pline[sl-1] == CONTINUE))
650 /* build one long line from several fragments - necessary for CMAP */
651 while (continuing(line))
653 status = cpp_read_line(&handle, STRLEN, line);
654 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
656 /* Since we depend on the '\' being present to continue to read, we copy line
657 * to a tmp string, strip the '\' from that string, and cat it to pline
659 tmp_line = strdup(line);
661 sl = strlen(tmp_line);
662 if ((sl > 0) && (tmp_line[sl-1] == CONTINUE))
664 tmp_line[sl-1] = ' ';
667 done = (status == eCPP_EOF);
670 if (status != eCPP_OK)
672 gmx_fatal(FARGS, cpp_error(&handle, status));
676 fprintf(out, "%s\n", line);
680 srenew(pline, strlen(pline)+strlen(tmp_line)+1);
681 strcat(pline, tmp_line);
685 /* skip trailing and leading spaces and comment text */
686 strip_comment (pline);
689 /* if there is something left... */
690 if ((int)strlen(pline) > 0)
692 if (pline[0] == OPENDIR)
694 /* A directive on this line: copy the directive
695 * without the brackets into dirstr, then
696 * skip spaces and tabs on either side of directive
698 dirstr = strdup((pline+1));
699 if ((dummy2 = strchr (dirstr, CLOSEDIR)) != NULL)
705 if ((newd = str2dir(dirstr)) == d_invalid)
707 sprintf(errbuf, "Invalid directive %s", dirstr);
708 warning_error(wi, errbuf);
712 /* Directive found */
715 fprintf(debug, "found directive '%s'\n", dir2str(newd));
717 if (DS_Check_Order (DS, newd))
724 /* we should print here which directives should have
725 been present, and which actually are */
726 gmx_fatal(FARGS, "%s\nInvalid order for directive %s",
727 cpp_error(&handle, eCPP_SYNTAX), dir2str(newd));
733 else if (d != d_invalid)
735 /* Not a directive, just a plain string
736 * use a gigantic switch to decode,
737 * if there is a valid directive!
744 gmx_fatal(FARGS, "%s\nFound a second defaults directive.\n",
745 cpp_error(&handle, eCPP_SYNTAX));
747 bReadDefaults = TRUE;
748 nscan = sscanf(pline, "%s%s%s%lf%lf%lf",
749 nb_str, comb_str, genpairs, &fLJ, &fQQ, &fPOW);
760 get_nbparm(nb_str, comb_str, &nb_funct, &comb, wi);
761 *combination_rule = comb;
764 bGenPairs = (gmx_strncasecmp(genpairs, "Y", 1) == 0);
765 if (nb_funct != eNBF_LJ && bGenPairs)
767 gmx_fatal(FARGS, "Generating pair parameters is only supported with LJ non-bonded interactions");
783 nb_funct = ifunc_index(d_nonbond_params, nb_funct);
787 push_at(symtab, atype, batype, pline, nb_funct,
788 &nbparam, bGenPairs ? &pair : NULL, wi);
792 push_bt(d, plist, 2, NULL, batype, pline, wi);
794 case d_constrainttypes:
795 push_bt(d, plist, 2, NULL, batype, pline, wi);
800 push_nbt(d, pair, atype, pline, F_LJ14, wi);
804 push_bt(d, plist, 2, atype, NULL, pline, wi);
808 push_bt(d, plist, 3, NULL, batype, pline, wi);
810 case d_dihedraltypes:
811 /* Special routine that can read both 2 and 4 atom dihedral definitions. */
812 push_dihedraltype(d, plist, batype, pline, wi);
815 case d_nonbond_params:
816 push_nbt(d, nbparam, atype, pline, nb_funct, wi);
821 srenew(block,nblock);
822 srenew(blockinfo,nblock);
823 blk0=&(block[nblock-1]);
824 bi0=&(blockinfo[nblock-1]);
827 push_molt(symtab,bi0,pline);
831 case d_implicit_genborn_params:
832 push_gb_params(atype, pline, wi);
835 case d_implicit_surface_params:
836 gmx_fatal(FARGS, "Implicit surface directive not supported yet.");
840 push_cmaptype(d, plist, 5, atype, batype, pline, wi);
848 if (opts->couple_moltype != NULL &&
849 (opts->couple_lam0 == ecouplamNONE ||
850 opts->couple_lam0 == ecouplamQ ||
851 opts->couple_lam1 == ecouplamNONE ||
852 opts->couple_lam1 == ecouplamQ))
854 dcatt = add_atomtype_decoupled(symtab, atype,
855 &nbparam, bGenPairs ? &pair : NULL);
857 ntype = get_atomtype_ntypes(atype);
858 ncombs = (ntype*(ntype+1))/2;
859 generate_nbparams(comb, nb_funct, &(plist[nb_funct]), atype, wi);
860 ncopy = copy_nbparams(nbparam, nb_funct, &(plist[nb_funct]),
862 fprintf(stderr, "Generated %d of the %d non-bonded parameter combinations\n", ncombs-ncopy, ncombs);
863 free_nbparam(nbparam, ntype);
866 gen_pairs(&(plist[nb_funct]), &(plist[F_LJ14]), fudgeLJ, comb);
867 ncopy = copy_nbparams(pair, nb_funct, &(plist[F_LJ14]),
869 fprintf(stderr, "Generated %d of the %d 1-4 parameter combinations\n", ncombs-ncopy, ncombs);
870 free_nbparam(pair, ntype);
872 /* Copy GBSA parameters to atomtype array? */
877 push_molt(symtab, &nmol, molinfo, pline, wi);
878 srenew(block2, nmol);
879 block2[nmol-1].nr = 0;
880 mi0 = &((*molinfo)[nmol-1]);
884 push_atom(symtab, &(mi0->cgs), &(mi0->atoms), atype, pline, &lastcg, wi);
888 push_bond(d, plist, mi0->plist, &(mi0->atoms), atype, pline, FALSE,
889 bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
899 case d_position_restraints:
900 case d_angle_restraints:
901 case d_angle_restraints_z:
902 case d_distance_restraints:
903 case d_orientation_restraints:
904 case d_dihedral_restraints:
907 case d_water_polarization:
908 case d_thole_polarization:
909 push_bond(d, plist, mi0->plist, &(mi0->atoms), atype, pline, TRUE,
910 bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
913 push_cmap(d, plist, mi0->plist, &(mi0->atoms), atype, pline, wi);
917 push_vsitesn(d, mi0->plist, &(mi0->atoms), pline, wi);
921 if (!block2[nmol-1].nr)
923 init_block2(&(block2[nmol-1]), mi0->atoms.nr);
925 push_excl(pline, &(block2[nmol-1]));
929 title = put_symtab(symtab, pline);
936 push_mol(nmol, *molinfo, pline, &whichmol, &nrcopies, wi);
937 mi0 = &((*molinfo)[whichmol]);
938 srenew(molb, nmolb+1);
939 molb[nmolb].type = whichmol;
940 molb[nmolb].nmol = nrcopies;
943 bCouple = (opts->couple_moltype != NULL &&
944 (gmx_strcasecmp("system", opts->couple_moltype) == 0 ||
945 gmx_strcasecmp(*(mi0->name), opts->couple_moltype) == 0));
948 nmol_couple += nrcopies;
951 if (mi0->atoms.nr == 0)
953 gmx_fatal(FARGS, "Molecule type '%s' contains no atoms",
957 "Excluding %d bonded neighbours molecule type '%s'\n",
958 mi0->nrexcl, *mi0->name);
959 sum_q(&mi0->atoms, nrcopies, &qt, &qBt);
960 if (!mi0->bProcessed)
963 generate_excl(mi0->nrexcl,
968 merge_excl(&(mi0->excls), &(block2[whichmol]));
969 done_block2(&(block2[whichmol]));
970 make_shake(mi0->plist, &mi0->atoms, opts->nshake);
974 /* nnb contains information about first,2nd,3rd bonded neighbors.
975 * Use this to generate GB 1-2,1-3,1-4 interactions when necessary.
977 if (bGenborn == TRUE)
979 generate_gb_exclusion_interactions(mi0, atype, &nnb);
986 convert_moltype_couple(mi0, dcatt, *fudgeQQ,
987 opts->couple_lam0, opts->couple_lam1,
989 nb_funct, &(plist[nb_funct]));
991 stupid_fill_block(&mi0->mols, mi0->atoms.nr, TRUE);
992 mi0->bProcessed = TRUE;
997 fprintf (stderr, "case: %d\n", d);
1007 status = cpp_close_file(&handle);
1008 if (status != eCPP_OK)
1010 gmx_fatal(FARGS, cpp_error(&handle, status));
1015 gmx_fio_fclose(out);
1018 if (opts->couple_moltype)
1020 if (nmol_couple == 0)
1022 gmx_fatal(FARGS, "Did not find any molecules of type '%s' for coupling",
1023 opts->couple_moltype);
1025 fprintf(stderr, "Coupling %d copies of molecule type '%s'\n",
1026 nmol_couple, opts->couple_moltype);
1029 /* this is not very clean, but fixes core dump on empty system name */
1032 title = put_symtab(symtab, "");
1034 if (fabs(qt) > 1e-4)
1036 sprintf(warn_buf, "System has non-zero total charge: %.6f\n%s\n", qt, floating_point_arithmetic_tip);
1037 warning_note(wi, warn_buf);
1039 if (fabs(qBt) > 1e-4 && !gmx_within_tol(qBt, qt, 1e-6))
1041 sprintf(warn_buf, "State B has non-zero total charge: %.6f\n%s\n", qBt, floating_point_arithmetic_tip);
1042 warning_note(wi, warn_buf);
1045 for (i = 0; i < nmol; i++)
1047 done_block2(&(block2[i]));
1051 done_bond_atomtype(&batype);
1061 char **do_top(gmx_bool bVerbose,
1062 const char *topfile,
1063 const char *topppfile,
1068 int *combination_rule,
1069 double *repulsion_power,
1071 gpp_atomtype_t atype,
1073 t_molinfo **molinfo,
1076 gmx_molblock_t **molblock,
1080 /* Tmpfile might contain a long path */
1081 const char *tmpfile;
1086 tmpfile = topppfile;
1095 printf("processing topology...\n");
1097 title = read_topol(topfile, tmpfile, opts->define, opts->include,
1098 symtab, atype, nrmols, molinfo,
1099 plist, combination_rule, repulsion_power,
1100 opts, fudgeQQ, nmolblock, molblock,
1101 ir->efep != efepNO, bGenborn, bZero, wi);
1102 if ((*combination_rule != eCOMB_GEOMETRIC) &&
1103 (ir->vdwtype == evdwUSER))
1105 warning(wi, "Using sigma/epsilon based combination rules with"
1106 " user supplied potential function may produce unwanted"
1114 static void generate_qmexcl_moltype(gmx_moltype_t *molt, unsigned char *grpnr,
1117 /* This routine expects molt->ilist to be of size F_NRE and ordered. */
1119 /* generates the exclusions between the individual QM atoms, as
1120 * these interactions should be handled by the QM subroutines and
1121 * not by the gromacs routines
1124 i, j, l, k = 0, jmax, qm_max = 0, qm_nr = 0, nratoms = 0, link_nr = 0, link_max = 0;
1126 *qm_arr = NULL, *link_arr = NULL, a1, a2, a3, a4, ftype = 0;
1132 *bQMMM, *blink, bexcl;
1134 /* First we search and select the QM atoms in an qm_arr array that
1135 * we use to create the exclusions.
1137 * we take the possibility into account that a user has defined more
1138 * than one QM group:
1140 * for that we also need to do this an ugly work-about just in case
1141 * the QM group contains the entire system...
1143 jmax = ir->opts.ngQM;
1145 /* we first search for all the QM atoms and put them in an array
1147 for (j = 0; j < jmax; j++)
1149 for (i = 0; i < molt->atoms.nr; i++)
1151 if (qm_nr >= qm_max)
1154 srenew(qm_arr, qm_max);
1156 if ((grpnr ? grpnr[i] : 0) == j)
1158 qm_arr[qm_nr++] = i;
1162 /* bQMMM[..] is an array containin TRUE/FALSE for atoms that are
1163 * QM/not QM. We first set all elements to false. Afterwards we use
1164 * the qm_arr to change the elements corresponding to the QM atoms
1167 snew(bQMMM, molt->atoms.nr);
1168 for (i = 0; i < molt->atoms.nr; i++)
1172 for (i = 0; i < qm_nr; i++)
1174 bQMMM[qm_arr[i]] = TRUE;
1177 /* We remove all bonded interactions (i.e. bonds,
1178 * angles, dihedrals, 1-4's), involving the QM atoms. The way they
1179 * are removed is as follows: if the interaction invloves 2 atoms,
1180 * it is removed if both atoms are QMatoms. If it involves 3 atoms,
1181 * it is removed if at least two of the atoms are QM atoms, if the
1182 * interaction involves 4 atoms, it is removed if there are at least
1183 * 2 QM atoms. Since this routine is called once before any forces
1184 * are computed, the top->idef.il[N].iatom[] array (see idef.h) can
1185 * be rewritten at this poitn without any problem. 25-9-2002 */
1187 /* first check weter we already have CONNBONDS: */
1188 if (molt->ilist[F_CONNBONDS].nr != 0)
1190 fprintf(stderr, "nr. of CONNBONDS present already: %d\n",
1191 molt->ilist[F_CONNBONDS].nr/3);
1192 ftype = molt->ilist[F_CONNBONDS].iatoms[0];
1193 k = molt->ilist[F_CONNBONDS].nr;
1195 /* now we delete all bonded interactions, except the ones describing
1196 * a chemical bond. These are converted to CONNBONDS
1198 for (i = 0; i < F_LJ; i++)
1200 if (i == F_CONNBONDS)
1204 nratoms = interaction_function[i].nratoms;
1206 while (j < molt->ilist[i].nr)
1212 a1 = molt->ilist[i].iatoms[j+1];
1213 a2 = molt->ilist[i].iatoms[j+2];
1214 bexcl = (bQMMM[a1] && bQMMM[a2]);
1215 /* a bonded beteen two QM atoms will be copied to the
1216 * CONNBONDS list, for reasons mentioned above
1218 if (bexcl && i < F_ANGLES)
1220 srenew(molt->ilist[F_CONNBONDS].iatoms, k+3);
1221 molt->ilist[F_CONNBONDS].nr += 3;
1222 molt->ilist[F_CONNBONDS].iatoms[k++] = ftype;
1223 molt->ilist[F_CONNBONDS].iatoms[k++] = a1;
1224 molt->ilist[F_CONNBONDS].iatoms[k++] = a2;
1228 a1 = molt->ilist[i].iatoms[j+1];
1229 a2 = molt->ilist[i].iatoms[j+2];
1230 a3 = molt->ilist[i].iatoms[j+3];
1231 bexcl = ((bQMMM[a1] && bQMMM[a2]) ||
1232 (bQMMM[a1] && bQMMM[a3]) ||
1233 (bQMMM[a2] && bQMMM[a3]));
1236 a1 = molt->ilist[i].iatoms[j+1];
1237 a2 = molt->ilist[i].iatoms[j+2];
1238 a3 = molt->ilist[i].iatoms[j+3];
1239 a4 = molt->ilist[i].iatoms[j+4];
1240 bexcl = ((bQMMM[a1] && bQMMM[a2] && bQMMM[a3]) ||
1241 (bQMMM[a1] && bQMMM[a2] && bQMMM[a4]) ||
1242 (bQMMM[a1] && bQMMM[a3] && bQMMM[a4]) ||
1243 (bQMMM[a2] && bQMMM[a3] && bQMMM[a4]));
1246 gmx_fatal(FARGS, "no such bonded interactions with %d atoms\n", nratoms);
1250 /* since the interaction involves QM atoms, these should be
1251 * removed from the MM ilist
1253 molt->ilist[i].nr -= (nratoms+1);
1254 for (l = j; l < molt->ilist[i].nr; l++)
1256 molt->ilist[i].iatoms[l] = molt->ilist[i].iatoms[l+(nratoms+1)];
1261 j += nratoms+1; /* the +1 is for the functype */
1265 /* Now, we search for atoms bonded to a QM atom because we also want
1266 * to exclude their nonbonded interactions with the QM atoms. The
1267 * reason for this is that this interaction is accounted for in the
1268 * linkatoms interaction with the QMatoms and would be counted
1271 for (i = 0; i < F_NRE; i++)
1276 while (j < molt->ilist[i].nr)
1278 a1 = molt->ilist[i].iatoms[j+1];
1279 a2 = molt->ilist[i].iatoms[j+2];
1280 if ((bQMMM[a1] && !bQMMM[a2]) || (!bQMMM[a1] && bQMMM[a2]))
1282 if (link_nr >= link_max)
1285 srenew(link_arr, link_max);
1289 link_arr[link_nr++] = a2;
1293 link_arr[link_nr++] = a1;
1300 snew(blink, molt->atoms.nr);
1301 for (i = 0; i < molt->atoms.nr; i++)
1305 for (i = 0; i < link_nr; i++)
1307 blink[link_arr[i]] = TRUE;
1309 /* creating the exclusion block for the QM atoms. Each QM atom has
1310 * as excluded elements all the other QMatoms (and itself).
1312 qmexcl.nr = molt->atoms.nr;
1313 qmexcl.nra = qm_nr*(qm_nr+link_nr)+link_nr*qm_nr;
1314 snew(qmexcl.index, qmexcl.nr+1);
1315 snew(qmexcl.a, qmexcl.nra);
1317 for (i = 0; i < qmexcl.nr; i++)
1319 qmexcl.index[i] = j;
1322 for (k = 0; k < qm_nr; k++)
1324 qmexcl.a[k+j] = qm_arr[k];
1326 for (k = 0; k < link_nr; k++)
1328 qmexcl.a[qm_nr+k+j] = link_arr[k];
1330 j += (qm_nr+link_nr);
1334 for (k = 0; k < qm_nr; k++)
1336 qmexcl.a[k+j] = qm_arr[k];
1341 qmexcl.index[qmexcl.nr] = j;
1343 /* and merging with the exclusions already present in sys.
1346 init_block2(&qmexcl2, molt->atoms.nr);
1347 b_to_b2(&qmexcl, &qmexcl2);
1348 merge_excl(&(molt->excls), &qmexcl2);
1349 done_block2(&qmexcl2);
1351 /* Finally, we also need to get rid of the pair interactions of the
1352 * classical atom bonded to the boundary QM atoms with the QMatoms,
1353 * as this interaction is already accounted for by the QM, so also
1354 * here we run the risk of double counting! We proceed in a similar
1355 * way as we did above for the other bonded interactions: */
1356 for (i = F_LJ14; i < F_COUL14; i++)
1358 nratoms = interaction_function[i].nratoms;
1360 while (j < molt->ilist[i].nr)
1363 a1 = molt->ilist[i].iatoms[j+1];
1364 a2 = molt->ilist[i].iatoms[j+2];
1365 bexcl = ((bQMMM[a1] && bQMMM[a2]) ||
1366 (blink[a1] && bQMMM[a2]) ||
1367 (bQMMM[a1] && blink[a2]));
1370 /* since the interaction involves QM atoms, these should be
1371 * removed from the MM ilist
1373 molt->ilist[i].nr -= (nratoms+1);
1374 for (k = j; k < molt->ilist[i].nr; k++)
1376 molt->ilist[i].iatoms[k] = molt->ilist[i].iatoms[k+(nratoms+1)];
1381 j += nratoms+1; /* the +1 is for the functype */
1390 } /* generate_qmexcl */
1392 void generate_qmexcl(gmx_mtop_t *sys, t_inputrec *ir, warninp_t wi)
1394 /* This routine expects molt->molt[m].ilist to be of size F_NRE and ordered.
1397 unsigned char *grpnr;
1398 int mb, mol, nat_mol, i, nr_mol_with_qm_atoms = 0;
1399 gmx_molblock_t *molb;
1402 grpnr = sys->groups.grpnr[egcQMMM];
1404 for (mb = 0; mb < sys->nmolblock; mb++)
1406 molb = &sys->molblock[mb];
1407 nat_mol = sys->moltype[molb->type].atoms.nr;
1408 for (mol = 0; mol < molb->nmol; mol++)
1411 for (i = 0; i < nat_mol; i++)
1413 if ((grpnr ? grpnr[i] : 0) < ir->opts.ngQM)
1420 nr_mol_with_qm_atoms++;
1423 /* We need to split this molblock */
1426 /* Split the molblock at this molecule */
1428 srenew(sys->molblock, sys->nmolblock);
1429 for (i = sys->nmolblock-2; i >= mb; i--)
1431 sys->molblock[i+1] = sys->molblock[i];
1433 sys->molblock[mb ].nmol = mol;
1434 sys->molblock[mb+1].nmol -= mol;
1436 molb = &sys->molblock[mb];
1440 /* Split the molblock after this molecule */
1442 srenew(sys->molblock, sys->nmolblock);
1443 molb = &sys->molblock[mb];
1444 for (i = sys->nmolblock-2; i >= mb; i--)
1446 sys->molblock[i+1] = sys->molblock[i];
1448 sys->molblock[mb ].nmol = 1;
1449 sys->molblock[mb+1].nmol -= 1;
1452 /* Add a moltype for the QMMM molecule */
1454 srenew(sys->moltype, sys->nmoltype);
1455 /* Copy the moltype struct */
1456 sys->moltype[sys->nmoltype-1] = sys->moltype[molb->type];
1457 /* Copy the exclusions to a new array, since this is the only
1458 * thing that needs to be modified for QMMM.
1460 copy_blocka(&sys->moltype[molb->type ].excls,
1461 &sys->moltype[sys->nmoltype-1].excls);
1462 /* Set the molecule type for the QMMM molblock */
1463 molb->type = sys->nmoltype - 1;
1465 generate_qmexcl_moltype(&sys->moltype[molb->type], grpnr, ir);
1473 if (nr_mol_with_qm_atoms > 1)
1475 /* generate a warning is there are QM atoms in different
1476 * topologies. In this case it is not possible at this stage to
1477 * mutualy exclude the non-bonded interactions via the
1478 * exclusions (AFAIK). Instead, the user is advised to use the
1479 * energy group exclusions in the mdp file
1482 "\nThe QM subsystem is divided over multiple topologies. "
1483 "The mutual non-bonded interactions cannot be excluded. "
1484 "There are two ways to achieve this:\n\n"
1485 "1) merge the topologies, such that the atoms of the QM "
1486 "subsystem are all present in one single topology file. "
1487 "In this case this warning will dissappear\n\n"
1488 "2) exclude the non-bonded interactions explicitly via the "
1489 "energygrp-excl option in the mdp file. if this is the case "
1490 "this warning may be ignored"