2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, 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/fileio/gmxfio.h"
52 #include "gromacs/gmxlib/ifunc.h"
53 #include "gromacs/gmxlib/warninp.h"
54 #include "gromacs/gmxpreprocess/gmxcpp.h"
55 #include "gromacs/gmxpreprocess/gpp_bond_atomtype.h"
56 #include "gromacs/gmxpreprocess/gpp_nextnb.h"
57 #include "gromacs/gmxpreprocess/grompp-impl.h"
58 #include "gromacs/gmxpreprocess/topdirs.h"
59 #include "gromacs/gmxpreprocess/toppush.h"
60 #include "gromacs/gmxpreprocess/topshake.h"
61 #include "gromacs/gmxpreprocess/toputil.h"
62 #include "gromacs/gmxpreprocess/vsite_parm.h"
63 #include "gromacs/math/units.h"
64 #include "gromacs/math/utilities.h"
65 #include "gromacs/mdlib/genborn.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/topology/block.h"
69 #include "gromacs/topology/symtab.h"
70 #include "gromacs/topology/topology.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/smalloc.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
;
123 nnn
= static_cast<int>(std::sqrt(static_cast<double>(ntp
)));
124 GMX_ASSERT(nnn
* nnn
== ntp
, "Number of pairs of generated non-bonded parameters should be a perfect square");
126 nrfpA
= interaction_function
[F_LJ14
].nrfpA
;
127 nrfpB
= interaction_function
[F_LJ14
].nrfpB
;
130 if ((nrfp
!= nrfpA
) || (nrfpA
!= nrfpB
))
132 gmx_incons("Number of force parameters in gen_pairs wrong");
135 fprintf(stderr
, "Generating 1-4 interactions: fudge = %g\n", fudge
);
138 fprintf(debug
, "Fudge factor for 1-4 interactions: %g\n", fudge
);
139 fprintf(debug
, "Holy Cow! there are %d types\n", ntp
);
141 snew(pairs
->param
, pairs
->nr
);
142 for (i
= 0; (i
< ntp
); i
++)
145 pairs
->param
[i
].a
[0] = i
/ nnn
;
146 pairs
->param
[i
].a
[1] = i
% nnn
;
147 /* Copy normal and FEP parameters and multiply by fudge factor */
151 for (j
= 0; (j
< nrfp
); j
++)
153 /* If we are using sigma/epsilon values, only the epsilon values
154 * should be scaled, but not sigma.
155 * The sigma values have even indices 0,2, etc.
157 if ((comb
== eCOMB_ARITHMETIC
|| comb
== eCOMB_GEOM_SIG_EPS
) && (j
%2 == 0))
166 pairs
->param
[i
].c
[j
] = scaling
*nbs
->param
[i
].c
[j
];
167 pairs
->param
[i
].c
[nrfp
+j
] = scaling
*nbs
->param
[i
].c
[j
];
172 double check_mol(gmx_mtop_t
*mtop
, warninp_t wi
)
175 int i
, mb
, nmol
, ri
, pt
;
180 /* Check mass and charge */
183 for (mb
= 0; mb
< mtop
->nmoltype
; mb
++)
185 atoms
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].atoms
;
186 nmol
= mtop
->molblock
[mb
].nmol
;
187 for (i
= 0; (i
< atoms
->nr
); i
++)
189 q
+= nmol
*atoms
->atom
[i
].q
;
190 m
= atoms
->atom
[i
].m
;
191 mB
= atoms
->atom
[i
].mB
;
192 pt
= atoms
->atom
[i
].ptype
;
193 /* If the particle is an atom or a nucleus it must have a mass,
194 * else, if it is a shell, a vsite or a bondshell it can have mass zero
196 if (((m
<= 0.0) || (mB
<= 0.0)) && ((pt
== eptAtom
) || (pt
== eptNucleus
)))
198 ri
= atoms
->atom
[i
].resind
;
199 sprintf(buf
, "atom %s (Res %s-%d) has mass %g (state A) / %g (state B)\n",
200 *(atoms
->atomname
[i
]),
201 *(atoms
->resinfo
[ri
].name
),
202 atoms
->resinfo
[ri
].nr
,
204 warning_error(wi
, buf
);
207 if (((m
!= 0) || (mB
!= 0)) && (pt
== eptVSite
))
209 ri
= atoms
->atom
[i
].resind
;
210 sprintf(buf
, "virtual site %s (Res %s-%d) has non-zero mass %g (state A) / %g (state B)\n"
211 " Check your topology.\n",
212 *(atoms
->atomname
[i
]),
213 *(atoms
->resinfo
[ri
].name
),
214 atoms
->resinfo
[ri
].nr
,
216 warning_error(wi
, buf
);
217 /* The following statements make LINCS break! */
218 /* atoms->atom[i].m=0; */
225 static void sum_q(t_atoms
*atoms
, int n
, double *qt
, double *qBt
)
233 for (i
= 0; i
< atoms
->nr
; i
++)
235 qmolA
+= atoms
->atom
[i
].q
;
236 qmolB
+= atoms
->atom
[i
].qB
;
238 /* Unfortunately an absolute comparison,
239 * but this avoids unnecessary warnings and gmx-users mails.
241 if (fabs(qmolA
) >= 1e-6 || fabs(qmolB
) >= 1e-6)
248 static void get_nbparm(char *nb_str
, char *comb_str
, int *nb
, int *comb
,
252 char warn_buf
[STRLEN
];
255 for (i
= 1; (i
< eNBF_NR
); i
++)
257 if (gmx_strcasecmp(nb_str
, enbf_names
[i
]) == 0)
264 *nb
= strtol(nb_str
, NULL
, 10);
266 if ((*nb
< 1) || (*nb
>= eNBF_NR
))
268 sprintf(warn_buf
, "Invalid nonbond function selector '%s' using %s",
269 nb_str
, enbf_names
[1]);
270 warning_error(wi
, warn_buf
);
274 for (i
= 1; (i
< eCOMB_NR
); i
++)
276 if (gmx_strcasecmp(comb_str
, ecomb_names
[i
]) == 0)
283 *comb
= strtol(comb_str
, NULL
, 10);
285 if ((*comb
< 1) || (*comb
>= eCOMB_NR
))
287 sprintf(warn_buf
, "Invalid combination rule selector '%s' using %s",
288 comb_str
, ecomb_names
[1]);
289 warning_error(wi
, warn_buf
);
294 static char ** cpp_opts(const char *define
, const char *include
,
299 const char *cppadds
[2];
300 char **cppopts
= NULL
;
301 const char *option
[2] = { "-D", "-I" };
302 const char *nopt
[2] = { "define", "include" };
306 char warn_buf
[STRLEN
];
309 cppadds
[1] = include
;
310 for (n
= 0; (n
< 2); n
++)
317 while ((*ptr
!= '\0') && isspace(*ptr
))
322 while ((*rptr
!= '\0') && !isspace(*rptr
))
330 strncpy(buf
, ptr
, len
);
331 if (strstr(ptr
, option
[n
]) != ptr
)
333 set_warning_line(wi
, "mdp file", -1);
334 sprintf(warn_buf
, "Malformed %s option %s", nopt
[n
], buf
);
335 warning(wi
, warn_buf
);
339 srenew(cppopts
, ++ncppopts
);
340 cppopts
[ncppopts
-1] = gmx_strdup(buf
);
348 srenew(cppopts
, ++ncppopts
);
349 cppopts
[ncppopts
-1] = NULL
;
356 find_gb_bondlength(t_params
*plist
, int ai
, int aj
, real
*length
)
363 for (i
= 0; i
< F_NRE
&& !found
; i
++)
367 for (j
= 0; j
< plist
[i
].nr
; j
++)
369 a1
= plist
[i
].param
[j
].a
[0];
370 a2
= plist
[i
].param
[j
].a
[1];
372 if ( (a1
== ai
&& a2
== aj
) || (a1
== aj
&& a2
== ai
))
374 /* Equilibrium bond distance */
375 *length
= plist
[i
].param
[j
].c
[0];
388 find_gb_anglelength(t_params
*plist
, int ai
, int ak
, real
*length
)
390 int i
, j
, a1
, a2
, a3
;
393 int status
, status1
, status2
;
397 for (i
= 0; i
< F_NRE
&& !found
; i
++)
401 for (j
= 0; j
< plist
[i
].nr
; j
++)
403 a1
= plist
[i
].param
[j
].a
[0];
404 a2
= plist
[i
].param
[j
].a
[1];
405 a3
= plist
[i
].param
[j
].a
[2];
407 /* We dont care what the middle atom is, but use it below */
408 if ( (a1
== ai
&& a3
== ak
) || (a1
== ak
&& a3
== ai
) )
410 /* Equilibrium bond distance */
411 a123
= plist
[i
].param
[j
].c
[0];
412 /* Use middle atom to find reference distances r12 and r23 */
413 status1
= find_gb_bondlength(plist
, a1
, a2
, &r12
);
414 status2
= find_gb_bondlength(plist
, a2
, a3
, &r23
);
416 if (status1
== 0 && status2
== 0)
418 /* cosine theorem to get r13 */
419 *length
= std::sqrt(r12
*r12
+r23
*r23
-(2*r12
*r23
*cos(a123
/RAD2DEG
)));
432 generate_gb_exclusion_interactions(t_molinfo
*mi
, gpp_atomtype_t atype
, t_nextnb
*nnb
)
434 int j
, n
, ai
, aj
, ti
, tj
;
439 real radiusi
, radiusj
;
440 real gb_radiusi
, gb_radiusj
;
441 real param_c2
, param_c4
;
447 for (n
= 1; n
<= nnb
->nrex
; n
++)
462 /* Put all higher-order exclusions into 1,4 list so we dont miss them */
469 for (ai
= 0; ai
< nnb
->nr
; ai
++)
471 ti
= at
->atom
[ai
].type
;
472 radiusi
= get_atomtype_radius(ti
, atype
);
473 gb_radiusi
= get_atomtype_gb_radius(ti
, atype
);
475 for (j
= 0; j
< nnb
->nrexcl
[ai
][n
]; j
++)
477 aj
= nnb
->a
[ai
][n
][j
];
479 /* Only add the interactions once */
482 tj
= at
->atom
[aj
].type
;
483 radiusj
= get_atomtype_radius(tj
, atype
);
484 gb_radiusj
= get_atomtype_gb_radius(tj
, atype
);
486 /* There is an exclusion of type "ftype" between atoms ai and aj */
490 /* Reference distance, not used for 1-4 interactions */
494 if (find_gb_bondlength(plist
, ai
, aj
, &distance
) != 0)
496 gmx_fatal(FARGS
, "Cannot find bond length for atoms %d-%d", ai
, aj
);
500 if (find_gb_anglelength(plist
, ai
, aj
, &distance
) != 0)
502 gmx_fatal(FARGS
, "Cannot find length for atoms %d-%d involved in angle", ai
, aj
);
509 /* Assign GB parameters */
511 param
.c
[0] = radiusi
+radiusj
;
512 /* Reference distance distance */
513 param
.c
[1] = distance
;
514 /* Still parameter */
515 param
.c
[2] = param_c2
;
517 param
.c
[3] = gb_radiusi
+gb_radiusj
;
519 param
.c
[4] = param_c4
;
521 /* Add it to the parameter list */
522 add_param_to_list(&plist
[ftype
], ¶m
);
531 static void make_atoms_sys(int nmolb
, const gmx_molblock_t
*molb
,
532 const t_molinfo
*molinfo
,
536 const t_atoms
*mol_atoms
;
541 for (mb
= 0; mb
< nmolb
; mb
++)
543 mol_atoms
= &molinfo
[molb
[mb
].type
].atoms
;
545 srenew(atoms
->atom
, atoms
->nr
+ molb
[mb
].nmol
*mol_atoms
->nr
);
547 for (m
= 0; m
< molb
[mb
].nmol
; m
++)
549 for (a
= 0; a
< mol_atoms
->nr
; a
++)
551 atoms
->atom
[atoms
->nr
++] = mol_atoms
->atom
[a
];
558 static char **read_topol(const char *infile
, const char *outfile
,
559 const char *define
, const char *include
,
561 gpp_atomtype_t atype
,
564 t_molinfo
**intermolecular_interactions
,
566 int *combination_rule
,
571 gmx_molblock_t
**molblock
,
579 char *pline
= NULL
, **title
= NULL
;
580 char line
[STRLEN
], errbuf
[256], comb_str
[256], nb_str
[256];
582 char *dirstr
, *dummy2
;
583 int nrcopies
, nmol
, nmolb
= 0, nscan
, ncombs
, ncopy
;
584 double fLJ
, fQQ
, fPOW
;
585 gmx_molblock_t
*molb
= NULL
;
586 t_molinfo
*mi0
= NULL
;
589 t_nbparam
**nbparam
, **pair
;
591 real fudgeLJ
= -1; /* Multiplication factor to generate 1-4 from LJ */
592 gmx_bool bReadDefaults
, bReadMolType
, bGenPairs
, bWarn_copy_A_B
;
593 double qt
= 0, qBt
= 0; /* total charge */
594 t_bond_atomtype batype
;
596 int dcatt
= -1, nmol_couple
;
597 /* File handling variables */
600 char *tmp_line
= NULL
;
601 char warn_buf
[STRLEN
];
602 const char *floating_point_arithmetic_tip
=
603 "Total charge should normally be an integer. See\n"
604 "http://www.gromacs.org/Documentation/Floating_Point_Arithmetic\n"
605 "for discussion on how close it should be to an integer.\n";
606 /* We need to open the output file before opening the input file,
607 * because cpp_open_file can change the current working directory.
611 out
= gmx_fio_fopen(outfile
, "w");
618 /* open input file */
619 status
= cpp_open_file(infile
, &handle
, cpp_opts(define
, include
, wi
));
622 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
625 /* some local variables */
626 DS_Init(&DS
); /* directive stack */
627 nmol
= 0; /* no molecules yet... */
628 d
= d_invalid
; /* first thing should be a directive */
629 nbparam
= NULL
; /* The temporary non-bonded matrix */
630 pair
= NULL
; /* The temporary pair interaction matrix */
631 block2
= NULL
; /* the extra exclusions */
634 *reppow
= 12.0; /* Default value for repulsion power */
636 *intermolecular_interactions
= NULL
;
638 /* Init the number of CMAP torsion angles and grid spacing */
639 plist
[F_CMAP
].grid_spacing
= 0;
640 plist
[F_CMAP
].nc
= 0;
642 bWarn_copy_A_B
= bFEP
;
644 batype
= init_bond_atomtype();
645 /* parse the actual file */
646 bReadDefaults
= FALSE
;
648 bReadMolType
= FALSE
;
653 status
= cpp_read_line(&handle
, STRLEN
, line
);
654 done
= (status
== eCPP_EOF
);
657 if (status
!= eCPP_OK
)
659 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
663 fprintf(out
, "%s\n", line
);
666 set_warning_line(wi
, cpp_cur_file(&handle
), cpp_cur_linenr(&handle
));
668 pline
= gmx_strdup(line
);
670 /* Strip trailing '\' from pline, if it exists */
672 if ((sl
> 0) && (pline
[sl
-1] == CONTINUE
))
677 /* build one long line from several fragments - necessary for CMAP */
678 while (continuing(line
))
680 status
= cpp_read_line(&handle
, STRLEN
, line
);
681 set_warning_line(wi
, cpp_cur_file(&handle
), cpp_cur_linenr(&handle
));
683 /* Since we depend on the '\' being present to continue to read, we copy line
684 * to a tmp string, strip the '\' from that string, and cat it to pline
686 tmp_line
= gmx_strdup(line
);
688 sl
= strlen(tmp_line
);
689 if ((sl
> 0) && (tmp_line
[sl
-1] == CONTINUE
))
691 tmp_line
[sl
-1] = ' ';
694 done
= (status
== eCPP_EOF
);
697 if (status
!= eCPP_OK
)
699 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
703 fprintf(out
, "%s\n", line
);
707 srenew(pline
, strlen(pline
)+strlen(tmp_line
)+1);
708 strcat(pline
, tmp_line
);
712 /* skip trailing and leading spaces and comment text */
713 strip_comment (pline
);
716 /* if there is something left... */
717 if ((int)strlen(pline
) > 0)
719 if (pline
[0] == OPENDIR
)
721 /* A directive on this line: copy the directive
722 * without the brackets into dirstr, then
723 * skip spaces and tabs on either side of directive
725 dirstr
= gmx_strdup((pline
+1));
726 if ((dummy2
= strchr (dirstr
, CLOSEDIR
)) != NULL
)
732 if ((newd
= str2dir(dirstr
)) == d_invalid
)
734 sprintf(errbuf
, "Invalid directive %s", dirstr
);
735 warning_error(wi
, errbuf
);
739 /* Directive found */
742 fprintf(debug
, "found directive '%s'\n", dir2str(newd
));
744 if (DS_Check_Order (DS
, newd
))
751 /* we should print here which directives should have
752 been present, and which actually are */
753 gmx_fatal(FARGS
, "%s\nInvalid order for directive %s",
754 cpp_error(&handle
, eCPP_SYNTAX
), dir2str(newd
));
758 if (d
== d_intermolecular_interactions
)
760 if (*intermolecular_interactions
== NULL
)
762 /* We (mis)use the moleculetype processing
763 * to process the intermolecular interactions
764 * by making a "molecule" of the size of the system.
766 snew(*intermolecular_interactions
, 1);
767 init_molinfo(*intermolecular_interactions
);
768 mi0
= *intermolecular_interactions
;
769 make_atoms_sys(nmolb
, molb
, *molinfo
,
776 else if (d
!= d_invalid
)
778 /* Not a directive, just a plain string
779 * use a gigantic switch to decode,
780 * if there is a valid directive!
787 gmx_fatal(FARGS
, "%s\nFound a second defaults directive.\n",
788 cpp_error(&handle
, eCPP_SYNTAX
));
790 bReadDefaults
= TRUE
;
791 nscan
= sscanf(pline
, "%s%s%s%lf%lf%lf",
792 nb_str
, comb_str
, genpairs
, &fLJ
, &fQQ
, &fPOW
);
803 get_nbparm(nb_str
, comb_str
, &nb_funct
, combination_rule
, wi
);
806 bGenPairs
= (gmx_strncasecmp(genpairs
, "Y", 1) == 0);
807 if (nb_funct
!= eNBF_LJ
&& bGenPairs
)
809 gmx_fatal(FARGS
, "Generating pair parameters is only supported with LJ non-bonded interactions");
825 nb_funct
= ifunc_index(d_nonbond_params
, nb_funct
);
829 push_at(symtab
, atype
, batype
, pline
, nb_funct
,
830 &nbparam
, bGenPairs
? &pair
: NULL
, wi
);
834 push_bt(d
, plist
, 2, NULL
, batype
, pline
, wi
);
836 case d_constrainttypes
:
837 push_bt(d
, plist
, 2, NULL
, batype
, pline
, wi
);
842 push_nbt(d
, pair
, atype
, pline
, F_LJ14
, wi
);
846 push_bt(d
, plist
, 2, atype
, NULL
, pline
, wi
);
850 push_bt(d
, plist
, 3, NULL
, batype
, pline
, wi
);
852 case d_dihedraltypes
:
853 /* Special routine that can read both 2 and 4 atom dihedral definitions. */
854 push_dihedraltype(d
, plist
, batype
, pline
, wi
);
857 case d_nonbond_params
:
858 push_nbt(d
, nbparam
, atype
, pline
, nb_funct
, wi
);
863 srenew(block,nblock);
864 srenew(blockinfo,nblock);
865 blk0=&(block[nblock-1]);
866 bi0=&(blockinfo[nblock-1]);
869 push_molt(symtab,bi0,pline);
873 case d_implicit_genborn_params
:
874 push_gb_params(atype
, pline
, wi
);
877 case d_implicit_surface_params
:
878 gmx_fatal(FARGS
, "Implicit surface directive not supported yet.");
882 push_cmaptype(d
, plist
, 5, atype
, batype
, pline
, wi
);
890 if (opts
->couple_moltype
!= NULL
&&
891 (opts
->couple_lam0
== ecouplamNONE
||
892 opts
->couple_lam0
== ecouplamQ
||
893 opts
->couple_lam1
== ecouplamNONE
||
894 opts
->couple_lam1
== ecouplamQ
))
896 dcatt
= add_atomtype_decoupled(symtab
, atype
,
897 &nbparam
, bGenPairs
? &pair
: NULL
);
899 ntype
= get_atomtype_ntypes(atype
);
900 ncombs
= (ntype
*(ntype
+1))/2;
901 generate_nbparams(*combination_rule
, nb_funct
, &(plist
[nb_funct
]), atype
, wi
);
902 ncopy
= copy_nbparams(nbparam
, nb_funct
, &(plist
[nb_funct
]),
904 fprintf(stderr
, "Generated %d of the %d non-bonded parameter combinations\n", ncombs
-ncopy
, ncombs
);
905 free_nbparam(nbparam
, ntype
);
908 gen_pairs(&(plist
[nb_funct
]), &(plist
[F_LJ14
]), fudgeLJ
, *combination_rule
);
909 ncopy
= copy_nbparams(pair
, nb_funct
, &(plist
[F_LJ14
]),
911 fprintf(stderr
, "Generated %d of the %d 1-4 parameter combinations\n", ncombs
-ncopy
, ncombs
);
912 free_nbparam(pair
, ntype
);
914 /* Copy GBSA parameters to atomtype array? */
919 push_molt(symtab
, &nmol
, molinfo
, pline
, wi
);
920 srenew(block2
, nmol
);
921 block2
[nmol
-1].nr
= 0;
922 mi0
= &((*molinfo
)[nmol
-1]);
926 push_atom(symtab
, &(mi0
->cgs
), &(mi0
->atoms
), atype
, pline
, &lastcg
, wi
);
930 push_bond(d
, plist
, mi0
->plist
, &(mi0
->atoms
), atype
, pline
, FALSE
,
931 bGenPairs
, *fudgeQQ
, bZero
, &bWarn_copy_A_B
, wi
);
934 push_bond(d
, plist
, mi0
->plist
, &(mi0
->atoms
), atype
, pline
, FALSE
,
935 FALSE
, 1.0, bZero
, &bWarn_copy_A_B
, wi
);
945 case d_position_restraints
:
946 case d_angle_restraints
:
947 case d_angle_restraints_z
:
948 case d_distance_restraints
:
949 case d_orientation_restraints
:
950 case d_dihedral_restraints
:
953 case d_water_polarization
:
954 case d_thole_polarization
:
955 push_bond(d
, plist
, mi0
->plist
, &(mi0
->atoms
), atype
, pline
, TRUE
,
956 bGenPairs
, *fudgeQQ
, bZero
, &bWarn_copy_A_B
, wi
);
959 push_cmap(d
, plist
, mi0
->plist
, &(mi0
->atoms
), atype
, pline
, wi
);
963 push_vsitesn(d
, mi0
->plist
, &(mi0
->atoms
), pline
, wi
);
966 GMX_ASSERT(block2
, "block2 must always be allocated so exclusions can be processed");
967 if (!block2
[nmol
-1].nr
)
969 init_block2(&(block2
[nmol
-1]), mi0
->atoms
.nr
);
971 push_excl(pline
, &(block2
[nmol
-1]));
975 title
= put_symtab(symtab
, pline
);
982 push_mol(nmol
, *molinfo
, pline
, &whichmol
, &nrcopies
, wi
);
983 mi0
= &((*molinfo
)[whichmol
]);
984 srenew(molb
, nmolb
+1);
985 molb
[nmolb
].type
= whichmol
;
986 molb
[nmolb
].nmol
= nrcopies
;
989 bCouple
= (opts
->couple_moltype
!= NULL
&&
990 (gmx_strcasecmp("system", opts
->couple_moltype
) == 0 ||
991 gmx_strcasecmp(*(mi0
->name
), opts
->couple_moltype
) == 0));
994 nmol_couple
+= nrcopies
;
997 if (mi0
->atoms
.nr
== 0)
999 gmx_fatal(FARGS
, "Molecule type '%s' contains no atoms",
1003 "Excluding %d bonded neighbours molecule type '%s'\n",
1004 mi0
->nrexcl
, *mi0
->name
);
1005 sum_q(&mi0
->atoms
, nrcopies
, &qt
, &qBt
);
1006 if (!mi0
->bProcessed
)
1009 generate_excl(mi0
->nrexcl
,
1014 merge_excl(&(mi0
->excls
), &(block2
[whichmol
]));
1015 done_block2(&(block2
[whichmol
]));
1016 make_shake(mi0
->plist
, &mi0
->atoms
, opts
->nshake
);
1020 /* nnb contains information about first,2nd,3rd bonded neighbors.
1021 * Use this to generate GB 1-2,1-3,1-4 interactions when necessary.
1023 if (bGenborn
== TRUE
)
1025 generate_gb_exclusion_interactions(mi0
, atype
, &nnb
);
1032 convert_moltype_couple(mi0
, dcatt
, *fudgeQQ
,
1033 opts
->couple_lam0
, opts
->couple_lam1
,
1035 nb_funct
, &(plist
[nb_funct
]));
1037 stupid_fill_block(&mi0
->mols
, mi0
->atoms
.nr
, TRUE
);
1038 mi0
->bProcessed
= TRUE
;
1043 fprintf (stderr
, "case: %d\n", (int)d
);
1044 gmx_incons("unknown directive");
1053 status
= cpp_close_file(&handle
);
1054 if (status
!= eCPP_OK
)
1056 gmx_fatal(FARGS
, cpp_error(&handle
, status
));
1061 gmx_fio_fclose(out
);
1064 if (opts
->couple_moltype
)
1066 if (nmol_couple
== 0)
1068 gmx_fatal(FARGS
, "Did not find any molecules of type '%s' for coupling",
1069 opts
->couple_moltype
);
1071 fprintf(stderr
, "Coupling %d copies of molecule type '%s'\n",
1072 nmol_couple
, opts
->couple_moltype
);
1075 /* this is not very clean, but fixes core dump on empty system name */
1078 title
= put_symtab(symtab
, "");
1080 if (fabs(qt
) > 1e-4)
1082 sprintf(warn_buf
, "System has non-zero total charge: %.6f\n%s\n", qt
, floating_point_arithmetic_tip
);
1083 warning_note(wi
, warn_buf
);
1085 if (fabs(qBt
) > 1e-4 && !gmx_within_tol(qBt
, qt
, 1e-6))
1087 sprintf(warn_buf
, "State B has non-zero total charge: %.6f\n%s\n", qBt
, floating_point_arithmetic_tip
);
1088 warning_note(wi
, warn_buf
);
1091 for (i
= 0; i
< nmol
; i
++)
1093 done_block2(&(block2
[i
]));
1097 done_bond_atomtype(&batype
);
1099 if (*intermolecular_interactions
!= NULL
)
1101 sfree(mi0
->atoms
.atom
);
1112 char **do_top(gmx_bool bVerbose
,
1113 const char *topfile
,
1114 const char *topppfile
,
1119 int *combination_rule
,
1120 double *repulsion_power
,
1122 gpp_atomtype_t atype
,
1124 t_molinfo
**molinfo
,
1125 t_molinfo
**intermolecular_interactions
,
1128 gmx_molblock_t
**molblock
,
1132 /* Tmpfile might contain a long path */
1133 const char *tmpfile
;
1138 tmpfile
= topppfile
;
1147 printf("processing topology...\n");
1149 title
= read_topol(topfile
, tmpfile
, opts
->define
, opts
->include
,
1151 nrmols
, molinfo
, intermolecular_interactions
,
1152 plist
, combination_rule
, repulsion_power
,
1153 opts
, fudgeQQ
, nmolblock
, molblock
,
1154 ir
->efep
!= efepNO
, bGenborn
, bZero
, wi
);
1155 if ((*combination_rule
!= eCOMB_GEOMETRIC
) &&
1156 (ir
->vdwtype
== evdwUSER
))
1158 warning(wi
, "Using sigma/epsilon based combination rules with"
1159 " user supplied potential function may produce unwanted"
1167 static void generate_qmexcl_moltype(gmx_moltype_t
*molt
, unsigned char *grpnr
,
1170 /* This routine expects molt->ilist to be of size F_NRE and ordered. */
1172 /* generates the exclusions between the individual QM atoms, as
1173 * these interactions should be handled by the QM subroutines and
1174 * not by the gromacs routines
1177 i
, j
, l
, k
= 0, jmax
, qm_max
= 0, qm_nr
= 0, nratoms
= 0, link_nr
= 0, link_max
= 0;
1179 *qm_arr
= NULL
, *link_arr
= NULL
, a1
, a2
, a3
, a4
, ftype
= 0;
1185 *bQMMM
, *blink
, bexcl
;
1187 /* First we search and select the QM atoms in an qm_arr array that
1188 * we use to create the exclusions.
1190 * we take the possibility into account that a user has defined more
1191 * than one QM group:
1193 * for that we also need to do this an ugly work-about just in case
1194 * the QM group contains the entire system...
1196 jmax
= ir
->opts
.ngQM
;
1198 /* we first search for all the QM atoms and put them in an array
1200 for (j
= 0; j
< jmax
; j
++)
1202 for (i
= 0; i
< molt
->atoms
.nr
; i
++)
1204 if (qm_nr
>= qm_max
)
1207 srenew(qm_arr
, qm_max
);
1209 if ((grpnr
? grpnr
[i
] : 0) == j
)
1211 qm_arr
[qm_nr
++] = i
;
1215 /* bQMMM[..] is an array containin TRUE/FALSE for atoms that are
1216 * QM/not QM. We first set all elements to false. Afterwards we use
1217 * the qm_arr to change the elements corresponding to the QM atoms
1220 snew(bQMMM
, molt
->atoms
.nr
);
1221 for (i
= 0; i
< molt
->atoms
.nr
; i
++)
1225 for (i
= 0; i
< qm_nr
; i
++)
1227 bQMMM
[qm_arr
[i
]] = TRUE
;
1230 /* We remove all bonded interactions (i.e. bonds,
1231 * angles, dihedrals, 1-4's), involving the QM atoms. The way they
1232 * are removed is as follows: if the interaction invloves 2 atoms,
1233 * it is removed if both atoms are QMatoms. If it involves 3 atoms,
1234 * it is removed if at least two of the atoms are QM atoms, if the
1235 * interaction involves 4 atoms, it is removed if there are at least
1236 * 2 QM atoms. Since this routine is called once before any forces
1237 * are computed, the top->idef.il[N].iatom[] array (see idef.h) can
1238 * be rewritten at this poitn without any problem. 25-9-2002 */
1240 /* first check weter we already have CONNBONDS: */
1241 if (molt
->ilist
[F_CONNBONDS
].nr
!= 0)
1243 fprintf(stderr
, "nr. of CONNBONDS present already: %d\n",
1244 molt
->ilist
[F_CONNBONDS
].nr
/3);
1245 ftype
= molt
->ilist
[F_CONNBONDS
].iatoms
[0];
1246 k
= molt
->ilist
[F_CONNBONDS
].nr
;
1248 /* now we delete all bonded interactions, except the ones describing
1249 * a chemical bond. These are converted to CONNBONDS
1251 for (i
= 0; i
< F_LJ
; i
++)
1253 if (i
== F_CONNBONDS
)
1257 nratoms
= interaction_function
[i
].nratoms
;
1259 while (j
< molt
->ilist
[i
].nr
)
1264 a1
= molt
->ilist
[i
].iatoms
[j
+1];
1265 a2
= molt
->ilist
[i
].iatoms
[j
+2];
1266 bexcl
= (bQMMM
[a1
] && bQMMM
[a2
]);
1267 /* a bonded beteen two QM atoms will be copied to the
1268 * CONNBONDS list, for reasons mentioned above
1270 if (bexcl
&& i
< F_ANGLES
)
1272 srenew(molt
->ilist
[F_CONNBONDS
].iatoms
, k
+3);
1273 molt
->ilist
[F_CONNBONDS
].nr
+= 3;
1274 molt
->ilist
[F_CONNBONDS
].iatoms
[k
++] = ftype
;
1275 molt
->ilist
[F_CONNBONDS
].iatoms
[k
++] = a1
;
1276 molt
->ilist
[F_CONNBONDS
].iatoms
[k
++] = a2
;
1280 a1
= molt
->ilist
[i
].iatoms
[j
+1];
1281 a2
= molt
->ilist
[i
].iatoms
[j
+2];
1282 a3
= molt
->ilist
[i
].iatoms
[j
+3];
1283 bexcl
= ((bQMMM
[a1
] && bQMMM
[a2
]) ||
1284 (bQMMM
[a1
] && bQMMM
[a3
]) ||
1285 (bQMMM
[a2
] && bQMMM
[a3
]));
1288 a1
= molt
->ilist
[i
].iatoms
[j
+1];
1289 a2
= molt
->ilist
[i
].iatoms
[j
+2];
1290 a3
= molt
->ilist
[i
].iatoms
[j
+3];
1291 a4
= molt
->ilist
[i
].iatoms
[j
+4];
1292 bexcl
= ((bQMMM
[a1
] && bQMMM
[a2
] && bQMMM
[a3
]) ||
1293 (bQMMM
[a1
] && bQMMM
[a2
] && bQMMM
[a4
]) ||
1294 (bQMMM
[a1
] && bQMMM
[a3
] && bQMMM
[a4
]) ||
1295 (bQMMM
[a2
] && bQMMM
[a3
] && bQMMM
[a4
]));
1298 gmx_fatal(FARGS
, "no such bonded interactions with %d atoms\n", nratoms
);
1303 /* since the interaction involves QM atoms, these should be
1304 * removed from the MM ilist
1306 molt
->ilist
[i
].nr
-= (nratoms
+1);
1307 for (l
= j
; l
< molt
->ilist
[i
].nr
; l
++)
1309 molt
->ilist
[i
].iatoms
[l
] = molt
->ilist
[i
].iatoms
[l
+(nratoms
+1)];
1314 j
+= nratoms
+1; /* the +1 is for the functype */
1318 /* Now, we search for atoms bonded to a QM atom because we also want
1319 * to exclude their nonbonded interactions with the QM atoms. The
1320 * reason for this is that this interaction is accounted for in the
1321 * linkatoms interaction with the QMatoms and would be counted
1324 for (i
= 0; i
< F_NRE
; i
++)
1329 while (j
< molt
->ilist
[i
].nr
)
1331 a1
= molt
->ilist
[i
].iatoms
[j
+1];
1332 a2
= molt
->ilist
[i
].iatoms
[j
+2];
1333 if ((bQMMM
[a1
] && !bQMMM
[a2
]) || (!bQMMM
[a1
] && bQMMM
[a2
]))
1335 if (link_nr
>= link_max
)
1338 srenew(link_arr
, link_max
);
1342 link_arr
[link_nr
++] = a2
;
1346 link_arr
[link_nr
++] = a1
;
1353 snew(blink
, molt
->atoms
.nr
);
1354 for (i
= 0; i
< molt
->atoms
.nr
; i
++)
1358 for (i
= 0; i
< link_nr
; i
++)
1360 blink
[link_arr
[i
]] = TRUE
;
1362 /* creating the exclusion block for the QM atoms. Each QM atom has
1363 * as excluded elements all the other QMatoms (and itself).
1365 qmexcl
.nr
= molt
->atoms
.nr
;
1366 qmexcl
.nra
= qm_nr
*(qm_nr
+link_nr
)+link_nr
*qm_nr
;
1367 snew(qmexcl
.index
, qmexcl
.nr
+1);
1368 snew(qmexcl
.a
, qmexcl
.nra
);
1370 for (i
= 0; i
< qmexcl
.nr
; i
++)
1372 qmexcl
.index
[i
] = j
;
1375 for (k
= 0; k
< qm_nr
; k
++)
1377 qmexcl
.a
[k
+j
] = qm_arr
[k
];
1379 for (k
= 0; k
< link_nr
; k
++)
1381 qmexcl
.a
[qm_nr
+k
+j
] = link_arr
[k
];
1383 j
+= (qm_nr
+link_nr
);
1387 for (k
= 0; k
< qm_nr
; k
++)
1389 qmexcl
.a
[k
+j
] = qm_arr
[k
];
1394 qmexcl
.index
[qmexcl
.nr
] = j
;
1396 /* and merging with the exclusions already present in sys.
1399 init_block2(&qmexcl2
, molt
->atoms
.nr
);
1400 b_to_b2(&qmexcl
, &qmexcl2
);
1401 merge_excl(&(molt
->excls
), &qmexcl2
);
1402 done_block2(&qmexcl2
);
1404 /* Finally, we also need to get rid of the pair interactions of the
1405 * classical atom bonded to the boundary QM atoms with the QMatoms,
1406 * as this interaction is already accounted for by the QM, so also
1407 * here we run the risk of double counting! We proceed in a similar
1408 * way as we did above for the other bonded interactions: */
1409 for (i
= F_LJ14
; i
< F_COUL14
; i
++)
1411 nratoms
= interaction_function
[i
].nratoms
;
1413 while (j
< molt
->ilist
[i
].nr
)
1415 a1
= molt
->ilist
[i
].iatoms
[j
+1];
1416 a2
= molt
->ilist
[i
].iatoms
[j
+2];
1417 bexcl
= ((bQMMM
[a1
] && bQMMM
[a2
]) ||
1418 (blink
[a1
] && bQMMM
[a2
]) ||
1419 (bQMMM
[a1
] && blink
[a2
]));
1422 /* since the interaction involves QM atoms, these should be
1423 * removed from the MM ilist
1425 molt
->ilist
[i
].nr
-= (nratoms
+1);
1426 for (k
= j
; k
< molt
->ilist
[i
].nr
; k
++)
1428 molt
->ilist
[i
].iatoms
[k
] = molt
->ilist
[i
].iatoms
[k
+(nratoms
+1)];
1433 j
+= nratoms
+1; /* the +1 is for the functype */
1442 } /* generate_qmexcl */
1444 void generate_qmexcl(gmx_mtop_t
*sys
, t_inputrec
*ir
, warninp_t wi
)
1446 /* This routine expects molt->molt[m].ilist to be of size F_NRE and ordered.
1449 unsigned char *grpnr
;
1450 int mb
, mol
, nat_mol
, i
, nr_mol_with_qm_atoms
= 0;
1451 gmx_molblock_t
*molb
;
1454 grpnr
= sys
->groups
.grpnr
[egcQMMM
];
1456 for (mb
= 0; mb
< sys
->nmolblock
; mb
++)
1458 molb
= &sys
->molblock
[mb
];
1459 nat_mol
= sys
->moltype
[molb
->type
].atoms
.nr
;
1460 for (mol
= 0; mol
< molb
->nmol
; mol
++)
1463 for (i
= 0; i
< nat_mol
; i
++)
1465 if ((grpnr
? grpnr
[i
] : 0) < ir
->opts
.ngQM
)
1472 nr_mol_with_qm_atoms
++;
1475 /* We need to split this molblock */
1478 /* Split the molblock at this molecule */
1480 srenew(sys
->molblock
, sys
->nmolblock
);
1481 for (i
= sys
->nmolblock
-2; i
>= mb
; i
--)
1483 sys
->molblock
[i
+1] = sys
->molblock
[i
];
1485 sys
->molblock
[mb
].nmol
= mol
;
1486 sys
->molblock
[mb
+1].nmol
-= mol
;
1488 molb
= &sys
->molblock
[mb
];
1492 /* Split the molblock after this molecule */
1494 srenew(sys
->molblock
, sys
->nmolblock
);
1495 molb
= &sys
->molblock
[mb
];
1496 for (i
= sys
->nmolblock
-2; i
>= mb
; i
--)
1498 sys
->molblock
[i
+1] = sys
->molblock
[i
];
1500 sys
->molblock
[mb
].nmol
= 1;
1501 sys
->molblock
[mb
+1].nmol
-= 1;
1504 /* Add a moltype for the QMMM molecule */
1506 srenew(sys
->moltype
, sys
->nmoltype
);
1507 /* Copy the moltype struct */
1508 sys
->moltype
[sys
->nmoltype
-1] = sys
->moltype
[molb
->type
];
1509 /* Copy the exclusions to a new array, since this is the only
1510 * thing that needs to be modified for QMMM.
1512 copy_blocka(&sys
->moltype
[molb
->type
].excls
,
1513 &sys
->moltype
[sys
->nmoltype
-1].excls
);
1514 /* Set the molecule type for the QMMM molblock */
1515 molb
->type
= sys
->nmoltype
- 1;
1517 generate_qmexcl_moltype(&sys
->moltype
[molb
->type
], grpnr
, ir
);
1525 if (nr_mol_with_qm_atoms
> 1)
1527 /* generate a warning is there are QM atoms in different
1528 * topologies. In this case it is not possible at this stage to
1529 * mutualy exclude the non-bonded interactions via the
1530 * exclusions (AFAIK). Instead, the user is advised to use the
1531 * energy group exclusions in the mdp file
1534 "\nThe QM subsystem is divided over multiple topologies. "
1535 "The mutual non-bonded interactions cannot be excluded. "
1536 "There are two ways to achieve this:\n\n"
1537 "1) merge the topologies, such that the atoms of the QM "
1538 "subsystem are all present in one single topology file. "
1539 "In this case this warning will dissappear\n\n"
1540 "2) exclude the non-bonded interactions explicitly via the "
1541 "energygrp-excl option in the mdp file. if this is the case "
1542 "this warning may be ignored"