1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
41 #include <sys/types.h>
60 #include "gmx_fatal.h"
62 #include "vsite_parm.h"
68 #include "gpp_nextnb.h"
72 #include "gpp_bond_atomtype.h"
76 #define CPPMARK '#' /* mark from cpp */
77 #define OPENDIR '[' /* starting sign for directive */
78 #define CLOSEDIR ']' /* ending sign for directive */
81 printf("line: %d, maxavail: %d\n", __LINE__, maxavail()); \
84 static void free_nbparam(t_nbparam **param, int nr)
88 for (i = 0; i < nr; i++)
95 static int copy_nbparams(t_nbparam **param, int ftype, t_params *plist, int nr)
103 for (i = 0; i < nr; i++)
105 for (j = 0; j <= i; j++)
107 if (param[i][j].bSet)
109 for (f = 0; f < nrfp; f++)
111 plist->param[nr*i+j].c[f] = param[i][j].c[f];
112 plist->param[nr*j+i].c[f] = param[i][j].c[f];
122 static void gen_pairs(t_params *nbs, t_params *pairs, real fudge, int comb, gmx_bool bVerbose)
124 int i, j, ntp, nrfp, nrfpA, nrfpB, nnn;
129 nrfpA = interaction_function[F_LJ14].nrfpA;
130 nrfpB = interaction_function[F_LJ14].nrfpB;
133 if ((nrfp != nrfpA) || (nrfpA != nrfpB))
135 gmx_incons("Number of force parameters in gen_pairs wrong");
138 fprintf(stderr, "Generating 1-4 interactions: fudge = %g\n", fudge);
141 fprintf(debug, "Fudge factor for 1-4 interactions: %g\n", fudge);
142 fprintf(debug, "Holy Cow! there are %d types\n", ntp);
144 snew(pairs->param, pairs->nr);
145 for (i = 0; (i < ntp); i++)
148 pairs->param[i].a[0] = i / nnn;
149 pairs->param[i].a[1] = i % nnn;
150 /* Copy normal and FEP parameters and multiply by fudge factor */
154 for (j = 0; (j < nrfp); j++)
156 /* If we are using sigma/epsilon values, only the epsilon values
157 * should be scaled, but not sigma.
158 * The sigma values have even indices 0,2, etc.
160 if ((comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) && (j%2 == 0))
169 pairs->param[i].c[j] = scaling*nbs->param[i].c[j];
170 pairs->param[i].c[nrfp+j] = scaling*nbs->param[i].c[j];
175 double check_mol(gmx_mtop_t *mtop, warninp_t wi)
178 int i, mb, nmol, ri, pt;
183 /* Check mass and charge */
186 for (mb = 0; mb < mtop->nmoltype; mb++)
188 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
189 nmol = mtop->molblock[mb].nmol;
190 for (i = 0; (i < atoms->nr); i++)
192 q += nmol*atoms->atom[i].q;
193 m = atoms->atom[i].m;
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) && ((pt == eptAtom) || (pt == eptNucleus)))
200 ri = atoms->atom[i].resind;
201 sprintf(buf, "atom %s (Res %s-%d) has mass %g\n",
202 *(atoms->atomname[i]),
203 *(atoms->resinfo[ri].name),
204 atoms->resinfo[ri].nr,
206 warning_error(wi, buf);
209 if ((m != 0) && (pt == eptVSite))
211 ri = atoms->atom[i].resind;
212 sprintf(buf, "virtual site %s (Res %s-%d) has non-zero mass %g\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,
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,
557 int i, sl, nb_funct, comb;
558 char *pline = NULL, **title = NULL;
559 char line[STRLEN], errbuf[256], comb_str[256], nb_str[256];
561 char *dirstr, *dummy2;
562 int nrcopies, nmol, nmolb = 0, nscan, ncombs, ncopy;
563 double fLJ, fQQ, fPOW;
564 gmx_molblock_t *molb = NULL;
565 t_topology *block = NULL;
566 t_molinfo *mi0 = NULL;
569 t_nbparam **nbparam, **pair;
571 real fudgeLJ = -1; /* Multiplication factor to generate 1-4 from LJ */
572 gmx_bool bReadDefaults, bReadMolType, bGenPairs, bWarn_copy_A_B;
573 double qt = 0, qBt = 0; /* total charge */
574 t_bond_atomtype batype;
576 int dcatt = -1, nmol_couple;
577 /* File handling variables */
580 char *tmp_line = NULL;
581 char warn_buf[STRLEN];
582 const char *floating_point_arithmetic_tip =
583 "Total charge should normally be an integer. See\n"
584 "http://www.gromacs.org/Documentation/Floating_Point_Arithmetic\n"
585 "for discussion on how close it should be to an integer.\n";
586 /* We need to open the output file before opening the input file,
587 * because cpp_open_file can change the current working directory.
591 out = gmx_fio_fopen(outfile, "w");
598 /* open input file */
599 status = cpp_open_file(infile, &handle, cpp_opts(define, include, infile, wi));
602 gmx_fatal(FARGS, cpp_error(&handle, status));
605 /* some local variables */
606 DS_Init(&DS); /* directive stack */
607 nmol = 0; /* no molecules yet... */
608 d = d_invalid; /* first thing should be a directive */
609 nbparam = NULL; /* The temporary non-bonded matrix */
610 pair = NULL; /* The temporary pair interaction matrix */
611 block2 = NULL; /* the extra exclusions */
613 *reppow = 12.0; /* Default value for repulsion power */
617 /* Init the number of CMAP torsion angles and grid spacing */
618 plist->grid_spacing = 0;
621 bWarn_copy_A_B = bFEP;
623 batype = init_bond_atomtype();
624 /* parse the actual file */
625 bReadDefaults = FALSE;
627 bReadMolType = FALSE;
632 status = cpp_read_line(&handle, STRLEN, line);
633 done = (status == eCPP_EOF);
636 if (status != eCPP_OK)
638 gmx_fatal(FARGS, cpp_error(&handle, status));
642 fprintf(out, "%s\n", line);
645 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
647 pline = strdup(line);
649 /* Strip trailing '\' from pline, if it exists */
651 if ((sl > 0) && (pline[sl-1] == CONTINUE))
656 /* build one long line from several fragments - necessary for CMAP */
657 while (continuing(line))
659 status = cpp_read_line(&handle, STRLEN, line);
660 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
662 /* Since we depend on the '\' being present to continue to read, we copy line
663 * to a tmp string, strip the '\' from that string, and cat it to pline
665 tmp_line = strdup(line);
667 sl = strlen(tmp_line);
668 if ((sl > 0) && (tmp_line[sl-1] == CONTINUE))
670 tmp_line[sl-1] = ' ';
673 done = (status == eCPP_EOF);
676 if (status != eCPP_OK)
678 gmx_fatal(FARGS, cpp_error(&handle, status));
682 fprintf(out, "%s\n", line);
686 srenew(pline, strlen(pline)+strlen(tmp_line)+1);
687 strcat(pline, tmp_line);
691 /* skip trailing and leading spaces and comment text */
692 strip_comment (pline);
695 /* if there is something left... */
696 if ((int)strlen(pline) > 0)
698 if (pline[0] == OPENDIR)
700 /* A directive on this line: copy the directive
701 * without the brackets into dirstr, then
702 * skip spaces and tabs on either side of directive
704 dirstr = strdup((pline+1));
705 if ((dummy2 = strchr (dirstr, CLOSEDIR)) != NULL)
711 if ((newd = str2dir(dirstr)) == d_invalid)
713 sprintf(errbuf, "Invalid directive %s", dirstr);
714 warning_error(wi, errbuf);
718 /* Directive found */
721 fprintf(debug, "found directive '%s'\n", dir2str(newd));
723 if (DS_Check_Order (DS, newd))
730 /* we should print here which directives should have
731 been present, and which actually are */
732 gmx_fatal(FARGS, "%s\nInvalid order for directive %s",
733 cpp_error(&handle, eCPP_SYNTAX), dir2str(newd));
739 else if (d != d_invalid)
741 /* Not a directive, just a plain string
742 * use a gigantic switch to decode,
743 * if there is a valid directive!
750 gmx_fatal(FARGS, "%s\nFound a second defaults directive.\n",
751 cpp_error(&handle, eCPP_SYNTAX));
753 bReadDefaults = TRUE;
754 nscan = sscanf(pline, "%s%s%s%lf%lf%lf",
755 nb_str, comb_str, genpairs, &fLJ, &fQQ, &fPOW);
766 get_nbparm(nb_str, comb_str, &nb_funct, &comb, wi);
767 *combination_rule = comb;
770 bGenPairs = (gmx_strncasecmp(genpairs, "Y", 1) == 0);
771 if (nb_funct != eNBF_LJ && bGenPairs)
773 gmx_fatal(FARGS, "Generating pair parameters is only supported with LJ non-bonded interactions");
789 nb_funct = ifunc_index(d_nonbond_params, nb_funct);
793 push_at(symtab, atype, batype, pline, nb_funct,
794 &nbparam, bGenPairs ? &pair : NULL, wi);
798 push_bt(d, plist, 2, NULL, batype, pline, wi);
800 case d_constrainttypes:
801 push_bt(d, plist, 2, NULL, batype, pline, wi);
806 push_nbt(d, pair, atype, pline, F_LJ14, wi);
810 push_bt(d, plist, 2, atype, NULL, pline, wi);
814 push_bt(d, plist, 3, NULL, batype, pline, wi);
816 case d_dihedraltypes:
817 /* Special routine that can read both 2 and 4 atom dihedral definitions. */
818 push_dihedraltype(d, plist, batype, pline, wi);
821 case d_nonbond_params:
822 push_nbt(d, nbparam, atype, pline, nb_funct, wi);
827 srenew(block,nblock);
828 srenew(blockinfo,nblock);
829 blk0=&(block[nblock-1]);
830 bi0=&(blockinfo[nblock-1]);
833 push_molt(symtab,bi0,pline);
837 case d_implicit_genborn_params:
838 push_gb_params(atype, pline, wi);
841 case d_implicit_surface_params:
842 gmx_fatal(FARGS, "Implicit surface directive not supported yet.");
846 push_cmaptype(d, plist, 5, atype, batype, pline, wi);
854 if (opts->couple_moltype != NULL &&
855 (opts->couple_lam0 == ecouplamNONE ||
856 opts->couple_lam0 == ecouplamQ ||
857 opts->couple_lam1 == ecouplamNONE ||
858 opts->couple_lam1 == ecouplamQ))
860 dcatt = add_atomtype_decoupled(symtab, atype,
861 &nbparam, bGenPairs ? &pair : NULL);
863 ntype = get_atomtype_ntypes(atype);
864 ncombs = (ntype*(ntype+1))/2;
865 generate_nbparams(comb, nb_funct, &(plist[nb_funct]), atype, wi);
866 ncopy = copy_nbparams(nbparam, nb_funct, &(plist[nb_funct]),
868 fprintf(stderr, "Generated %d of the %d non-bonded parameter combinations\n", ncombs-ncopy, ncombs);
869 free_nbparam(nbparam, ntype);
872 gen_pairs(&(plist[nb_funct]), &(plist[F_LJ14]), fudgeLJ, comb, bVerbose);
873 ncopy = copy_nbparams(pair, nb_funct, &(plist[F_LJ14]),
875 fprintf(stderr, "Generated %d of the %d 1-4 parameter combinations\n", ncombs-ncopy, ncombs);
876 free_nbparam(pair, ntype);
878 /* Copy GBSA parameters to atomtype array? */
883 push_molt(symtab, &nmol, molinfo, pline, wi);
884 srenew(block2, nmol);
885 block2[nmol-1].nr = 0;
886 mi0 = &((*molinfo)[nmol-1]);
890 push_atom(symtab, &(mi0->cgs), &(mi0->atoms), atype, pline, &lastcg, wi);
894 push_bond(d, plist, mi0->plist, &(mi0->atoms), atype, pline, FALSE,
895 bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
905 case d_position_restraints:
906 case d_angle_restraints:
907 case d_angle_restraints_z:
908 case d_distance_restraints:
909 case d_orientation_restraints:
910 case d_dihedral_restraints:
913 case d_water_polarization:
914 case d_thole_polarization:
915 push_bond(d, plist, mi0->plist, &(mi0->atoms), atype, pline, TRUE,
916 bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
919 push_cmap(d, plist, mi0->plist, &(mi0->atoms), atype, pline,
920 &bWarn_copy_A_B, wi);
924 push_vsitesn(d, plist, mi0->plist, &(mi0->atoms), atype, pline, wi);
928 if (!block2[nmol-1].nr)
930 init_block2(&(block2[nmol-1]), mi0->atoms.nr);
932 push_excl(pline, &(block2[nmol-1]));
936 title = put_symtab(symtab, pline);
943 push_mol(nmol, *molinfo, pline, &whichmol, &nrcopies, wi);
944 mi0 = &((*molinfo)[whichmol]);
945 srenew(molb, nmolb+1);
946 molb[nmolb].type = whichmol;
947 molb[nmolb].nmol = nrcopies;
950 bCouple = (opts->couple_moltype != NULL &&
951 (gmx_strcasecmp("system", opts->couple_moltype) == 0 ||
952 gmx_strcasecmp(*(mi0->name), opts->couple_moltype) == 0));
955 nmol_couple += nrcopies;
958 if (mi0->atoms.nr == 0)
960 gmx_fatal(FARGS, "Molecule type '%s' contains no atoms",
964 "Excluding %d bonded neighbours molecule type '%s'\n",
965 mi0->nrexcl, *mi0->name);
966 sum_q(&mi0->atoms, nrcopies, &qt, &qBt);
967 if (!mi0->bProcessed)
970 generate_excl(mi0->nrexcl,
975 merge_excl(&(mi0->excls), &(block2[whichmol]));
976 done_block2(&(block2[whichmol]));
977 make_shake(mi0->plist, &mi0->atoms, atype, opts->nshake);
981 /* nnb contains information about first,2nd,3rd bonded neighbors.
982 * Use this to generate GB 1-2,1-3,1-4 interactions when necessary.
984 if (bGenborn == TRUE)
986 generate_gb_exclusion_interactions(mi0, atype, &nnb);
993 convert_moltype_couple(mi0, dcatt, *fudgeQQ,
994 opts->couple_lam0, opts->couple_lam1,
996 nb_funct, &(plist[nb_funct]));
998 stupid_fill_block(&mi0->mols, mi0->atoms.nr, TRUE);
999 mi0->bProcessed = TRUE;
1004 fprintf (stderr, "case: %d\n", d);
1014 status = cpp_close_file(&handle);
1015 if (status != eCPP_OK)
1017 gmx_fatal(FARGS, cpp_error(&handle, status));
1021 gmx_fio_fclose(out);
1024 if (opts->couple_moltype)
1026 if (nmol_couple == 0)
1028 gmx_fatal(FARGS, "Did not find any molecules of type '%s' for coupling",
1029 opts->couple_moltype);
1031 fprintf(stderr, "Coupling %d copies of molecule type '%s'\n",
1032 nmol_couple, opts->couple_moltype);
1035 /* this is not very clean, but fixes core dump on empty system name */
1038 title = put_symtab(symtab, "");
1040 if (fabs(qt) > 1e-4)
1042 sprintf(warn_buf, "System has non-zero total charge: %.6f\n%s\n", qt, floating_point_arithmetic_tip);
1043 warning_note(wi, warn_buf);
1045 if (fabs(qBt) > 1e-4 && !gmx_within_tol(qBt, qt, 1e-6))
1047 sprintf(warn_buf, "State B has non-zero total charge: %.6f\n%s\n", qBt, floating_point_arithmetic_tip);
1048 warning_note(wi, warn_buf);
1051 for (i = 0; i < nmol; i++)
1053 done_block2(&(block2[i]));
1057 done_bond_atomtype(&batype);
1067 char **do_top(gmx_bool bVerbose,
1068 const char *topfile,
1069 const char *topppfile,
1074 int *combination_rule,
1075 double *repulsion_power,
1077 gpp_atomtype_t atype,
1079 t_molinfo **molinfo,
1082 gmx_molblock_t **molblock,
1086 /* Tmpfile might contain a long path */
1087 const char *tmpfile;
1092 tmpfile = topppfile;
1101 printf("processing topology...\n");
1103 title = read_topol(topfile, tmpfile, opts->define, opts->include,
1104 symtab, atype, nrmols, molinfo,
1105 plist, combination_rule, repulsion_power,
1106 opts, fudgeQQ, nmolblock, molblock,
1107 ir->efep != efepNO, bGenborn, bZero, bVerbose,
1109 if ((*combination_rule != eCOMB_GEOMETRIC) &&
1110 (ir->vdwtype == evdwUSER))
1112 warning(wi, "Using sigma/epsilon based combination rules with"
1113 " user supplied potential function may produce unwanted"
1121 static void generate_qmexcl_moltype(gmx_moltype_t *molt, unsigned char *grpnr,
1124 /* This routine expects molt->ilist to be of size F_NRE and ordered. */
1126 /* generates the exclusions between the individual QM atoms, as
1127 * these interactions should be handled by the QM subroutines and
1128 * not by the gromacs routines
1131 i, j, l, k = 0, jmax, qm_max = 0, qm_nr = 0, nratoms = 0, link_nr = 0, link_max = 0;
1133 *qm_arr = NULL, *link_arr = NULL, a1, a2, a3, a4, ftype = 0;
1139 *bQMMM, *blink, bexcl;
1141 /* First we search and select the QM atoms in an qm_arr array that
1142 * we use to create the exclusions.
1144 * we take the possibility into account that a user has defined more
1145 * than one QM group:
1147 * for that we also need to do this an ugly work-about just in case
1148 * the QM group contains the entire system...
1150 jmax = ir->opts.ngQM;
1152 /* we first search for all the QM atoms and put them in an array
1154 for (j = 0; j < jmax; j++)
1156 for (i = 0; i < molt->atoms.nr; i++)
1158 if (qm_nr >= qm_max)
1161 srenew(qm_arr, qm_max);
1163 if ((grpnr ? grpnr[i] : 0) == j)
1165 qm_arr[qm_nr++] = i;
1169 /* bQMMM[..] is an array containin TRUE/FALSE for atoms that are
1170 * QM/not QM. We first set all elements to false. Afterwards we use
1171 * the qm_arr to change the elements corresponding to the QM atoms
1174 snew(bQMMM, molt->atoms.nr);
1175 for (i = 0; i < molt->atoms.nr; i++)
1179 for (i = 0; i < qm_nr; i++)
1181 bQMMM[qm_arr[i]] = TRUE;
1184 /* We remove all bonded interactions (i.e. bonds,
1185 * angles, dihedrals, 1-4's), involving the QM atoms. The way they
1186 * are removed is as follows: if the interaction invloves 2 atoms,
1187 * it is removed if both atoms are QMatoms. If it involves 3 atoms,
1188 * it is removed if at least two of the atoms are QM atoms, if the
1189 * interaction involves 4 atoms, it is removed if there are at least
1190 * 2 QM atoms. Since this routine is called once before any forces
1191 * are computed, the top->idef.il[N].iatom[] array (see idef.h) can
1192 * be rewritten at this poitn without any problem. 25-9-2002 */
1194 /* first check weter we already have CONNBONDS: */
1195 if (molt->ilist[F_CONNBONDS].nr != 0)
1197 fprintf(stderr, "nr. of CONNBONDS present already: %d\n",
1198 molt->ilist[F_CONNBONDS].nr/3);
1199 ftype = molt->ilist[F_CONNBONDS].iatoms[0];
1200 k = molt->ilist[F_CONNBONDS].nr;
1202 /* now we delete all bonded interactions, except the ones describing
1203 * a chemical bond. These are converted to CONNBONDS
1205 for (i = 0; i < F_LJ; i++)
1207 if (i == F_CONNBONDS)
1211 nratoms = interaction_function[i].nratoms;
1213 while (j < molt->ilist[i].nr)
1219 a1 = molt->ilist[i].iatoms[j+1];
1220 a2 = molt->ilist[i].iatoms[j+2];
1221 bexcl = (bQMMM[a1] && bQMMM[a2]);
1222 /* a bonded beteen two QM atoms will be copied to the
1223 * CONNBONDS list, for reasons mentioned above
1225 if (bexcl && i < F_ANGLES)
1227 srenew(molt->ilist[F_CONNBONDS].iatoms, k+3);
1228 molt->ilist[F_CONNBONDS].nr += 3;
1229 molt->ilist[F_CONNBONDS].iatoms[k++] = ftype;
1230 molt->ilist[F_CONNBONDS].iatoms[k++] = a1;
1231 molt->ilist[F_CONNBONDS].iatoms[k++] = a2;
1235 a1 = molt->ilist[i].iatoms[j+1];
1236 a2 = molt->ilist[i].iatoms[j+2];
1237 a3 = molt->ilist[i].iatoms[j+3];
1238 bexcl = ((bQMMM[a1] && bQMMM[a2]) ||
1239 (bQMMM[a1] && bQMMM[a3]) ||
1240 (bQMMM[a2] && bQMMM[a3]));
1243 a1 = molt->ilist[i].iatoms[j+1];
1244 a2 = molt->ilist[i].iatoms[j+2];
1245 a3 = molt->ilist[i].iatoms[j+3];
1246 a4 = molt->ilist[i].iatoms[j+4];
1247 bexcl = ((bQMMM[a1] && bQMMM[a2] && bQMMM[a3]) ||
1248 (bQMMM[a1] && bQMMM[a2] && bQMMM[a4]) ||
1249 (bQMMM[a1] && bQMMM[a3] && bQMMM[a4]) ||
1250 (bQMMM[a2] && bQMMM[a3] && bQMMM[a4]));
1253 gmx_fatal(FARGS, "no such bonded interactions with %d atoms\n", nratoms);
1257 /* since the interaction involves QM atoms, these should be
1258 * removed from the MM ilist
1260 molt->ilist[i].nr -= (nratoms+1);
1261 for (l = j; l < molt->ilist[i].nr; l++)
1263 molt->ilist[i].iatoms[l] = molt->ilist[i].iatoms[l+(nratoms+1)];
1268 j += nratoms+1; /* the +1 is for the functype */
1272 /* Now, we search for atoms bonded to a QM atom because we also want
1273 * to exclude their nonbonded interactions with the QM atoms. The
1274 * reason for this is that this interaction is accounted for in the
1275 * linkatoms interaction with the QMatoms and would be counted
1278 for (i = 0; i < F_NRE; i++)
1283 while (j < molt->ilist[i].nr)
1285 a1 = molt->ilist[i].iatoms[j+1];
1286 a2 = molt->ilist[i].iatoms[j+2];
1287 if ((bQMMM[a1] && !bQMMM[a2]) || (!bQMMM[a1] && bQMMM[a2]))
1289 if (link_nr >= link_max)
1292 srenew(link_arr, link_max);
1296 link_arr[link_nr++] = a2;
1300 link_arr[link_nr++] = a1;
1307 snew(blink, molt->atoms.nr);
1308 for (i = 0; i < molt->atoms.nr; i++)
1312 for (i = 0; i < link_nr; i++)
1314 blink[link_arr[i]] = TRUE;
1316 /* creating the exclusion block for the QM atoms. Each QM atom has
1317 * as excluded elements all the other QMatoms (and itself).
1319 qmexcl.nr = molt->atoms.nr;
1320 qmexcl.nra = qm_nr*(qm_nr+link_nr)+link_nr*qm_nr;
1321 snew(qmexcl.index, qmexcl.nr+1);
1322 snew(qmexcl.a, qmexcl.nra);
1324 for (i = 0; i < qmexcl.nr; i++)
1326 qmexcl.index[i] = j;
1329 for (k = 0; k < qm_nr; k++)
1331 qmexcl.a[k+j] = qm_arr[k];
1333 for (k = 0; k < link_nr; k++)
1335 qmexcl.a[qm_nr+k+j] = link_arr[k];
1337 j += (qm_nr+link_nr);
1341 for (k = 0; k < qm_nr; k++)
1343 qmexcl.a[k+j] = qm_arr[k];
1348 qmexcl.index[qmexcl.nr] = j;
1350 /* and merging with the exclusions already present in sys.
1353 init_block2(&qmexcl2, molt->atoms.nr);
1354 b_to_b2(&qmexcl, &qmexcl2);
1355 merge_excl(&(molt->excls), &qmexcl2);
1356 done_block2(&qmexcl2);
1358 /* Finally, we also need to get rid of the pair interactions of the
1359 * classical atom bonded to the boundary QM atoms with the QMatoms,
1360 * as this interaction is already accounted for by the QM, so also
1361 * here we run the risk of double counting! We proceed in a similar
1362 * way as we did above for the other bonded interactions: */
1363 for (i = F_LJ14; i < F_COUL14; i++)
1365 nratoms = interaction_function[i].nratoms;
1367 while (j < molt->ilist[i].nr)
1370 a1 = molt->ilist[i].iatoms[j+1];
1371 a2 = molt->ilist[i].iatoms[j+2];
1372 bexcl = ((bQMMM[a1] && bQMMM[a2]) ||
1373 (blink[a1] && bQMMM[a2]) ||
1374 (bQMMM[a1] && blink[a2]));
1377 /* since the interaction involves QM atoms, these should be
1378 * removed from the MM ilist
1380 molt->ilist[i].nr -= (nratoms+1);
1381 for (k = j; k < molt->ilist[i].nr; k++)
1383 molt->ilist[i].iatoms[k] = molt->ilist[i].iatoms[k+(nratoms+1)];
1388 j += nratoms+1; /* the +1 is for the functype */
1397 } /* generate_qmexcl */
1399 void generate_qmexcl(gmx_mtop_t *sys, t_inputrec *ir, warninp_t wi)
1401 /* This routine expects molt->molt[m].ilist to be of size F_NRE and ordered.
1404 unsigned char *grpnr;
1405 int mb, mol, nat_mol, i, nr_mol_with_qm_atoms = 0;
1406 gmx_molblock_t *molb;
1409 grpnr = sys->groups.grpnr[egcQMMM];
1411 for (mb = 0; mb < sys->nmolblock; mb++)
1413 molb = &sys->molblock[mb];
1414 nat_mol = sys->moltype[molb->type].atoms.nr;
1415 for (mol = 0; mol < molb->nmol; mol++)
1418 for (i = 0; i < nat_mol; i++)
1420 if ((grpnr ? grpnr[i] : 0) < ir->opts.ngQM)
1427 nr_mol_with_qm_atoms++;
1430 /* We need to split this molblock */
1433 /* Split the molblock at this molecule */
1435 srenew(sys->molblock, sys->nmolblock);
1436 for (i = sys->nmolblock-2; i >= mb; i--)
1438 sys->molblock[i+1] = sys->molblock[i];
1440 sys->molblock[mb ].nmol = mol;
1441 sys->molblock[mb+1].nmol -= mol;
1443 molb = &sys->molblock[mb];
1447 /* Split the molblock after this molecule */
1449 srenew(sys->molblock, sys->nmolblock);
1450 molb = &sys->molblock[mb];
1451 for (i = sys->nmolblock-2; i >= mb; i--)
1453 sys->molblock[i+1] = sys->molblock[i];
1455 sys->molblock[mb ].nmol = 1;
1456 sys->molblock[mb+1].nmol -= 1;
1459 /* Add a moltype for the QMMM molecule */
1461 srenew(sys->moltype, sys->nmoltype);
1462 /* Copy the moltype struct */
1463 sys->moltype[sys->nmoltype-1] = sys->moltype[molb->type];
1464 /* Copy the exclusions to a new array, since this is the only
1465 * thing that needs to be modified for QMMM.
1467 copy_blocka(&sys->moltype[molb->type ].excls,
1468 &sys->moltype[sys->nmoltype-1].excls);
1469 /* Set the molecule type for the QMMM molblock */
1470 molb->type = sys->nmoltype - 1;
1472 generate_qmexcl_moltype(&sys->moltype[molb->type], grpnr, ir);
1480 if (nr_mol_with_qm_atoms > 1)
1482 /* generate a warning is there are QM atoms in different
1483 * topologies. In this case it is not possible at this stage to
1484 * mutualy exclude the non-bonded interactions via the
1485 * exclusions (AFAIK). Instead, the user is advised to use the
1486 * energy group exclusions in the mdp file
1489 "\nThe QM subsystem is divided over multiple topologies. "
1490 "The mutual non-bonded interactions cannot be excluded. "
1491 "There are two ways to achieve this:\n\n"
1492 "1) merge the topologies, such that the atoms of the QM "
1493 "subsystem are all present in one single topology file. "
1494 "In this case this warning will dissappear\n\n"
1495 "2) exclude the non-bonded interactions explicitly via the "
1496 "energygrp-excl option in the mdp file. if this is the case "
1497 "this warning may be ignored"