Merge branch 'release-4-6' into release-5-0
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / topio.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include <sys/types.h>
43 #include <stdio.h>
44 #include <string.h>
45 #include <errno.h>
46 #include <ctype.h>
47 #include <assert.h>
48
49 #include "gromacs/fileio/futil.h"
50 #include "sysstuff.h"
51 #include "typedefs.h"
52 #include "gromacs/utility/smalloc.h"
53 #include "macros.h"
54 #include "gromacs/fileio/gmxfio.h"
55 #include "txtdump.h"
56 #include "physics.h"
57 #include "macros.h"
58 #include "names.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "symtab.h"
61 #include "gmx_fatal.h"
62 #include "warninp.h"
63 #include "vsite_parm.h"
64
65 #include "grompp-impl.h"
66 #include "toputil.h"
67 #include "toppush.h"
68 #include "topdirs.h"
69 #include "gpp_nextnb.h"
70 #include "topio.h"
71 #include "topshake.h"
72 #include "gmxcpp.h"
73 #include "gpp_bond_atomtype.h"
74 #include "genborn.h"
75 #include "gromacs/math/utilities.h"
76
77 #define OPENDIR     '[' /* starting sign for directive */
78 #define CLOSEDIR    ']' /* ending sign for directive   */
79
80 static void free_nbparam(t_nbparam **param, int nr)
81 {
82     int i;
83
84     for (i = 0; i < nr; i++)
85     {
86         sfree(param[i]);
87     }
88     sfree(param);
89 }
90
91 static int copy_nbparams(t_nbparam **param, int ftype, t_params *plist, int nr)
92 {
93     int i, j, f;
94     int nrfp, ncopy;
95
96     nrfp = NRFP(ftype);
97
98     ncopy = 0;
99     for (i = 0; i < nr; i++)
100     {
101         for (j = 0; j <= i; j++)
102         {
103             if (param[i][j].bSet)
104             {
105                 for (f = 0; f < nrfp; f++)
106                 {
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];
109                 }
110                 ncopy++;
111             }
112         }
113     }
114
115     return ncopy;
116 }
117
118 static void gen_pairs(t_params *nbs, t_params *pairs, real fudge, int comb)
119 {
120     int     i, j, ntp, nrfp, nrfpA, nrfpB, nnn;
121     real    scaling;
122     ntp       = nbs->nr;
123     nnn       = sqrt(ntp);
124     nrfp      = NRFP(F_LJ);
125     nrfpA     = interaction_function[F_LJ14].nrfpA;
126     nrfpB     = interaction_function[F_LJ14].nrfpB;
127     pairs->nr = ntp;
128
129     if ((nrfp  != nrfpA) || (nrfpA != nrfpB))
130     {
131         gmx_incons("Number of force parameters in gen_pairs wrong");
132     }
133
134     fprintf(stderr, "Generating 1-4 interactions: fudge = %g\n", fudge);
135     if (debug)
136     {
137         fprintf(debug, "Fudge factor for 1-4 interactions: %g\n", fudge);
138         fprintf(debug, "Holy Cow! there are %d types\n", ntp);
139     }
140     snew(pairs->param, pairs->nr);
141     for (i = 0; (i < ntp); i++)
142     {
143         /* Copy param.a */
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 */
147
148
149
150         for (j = 0; (j < nrfp); j++)
151         {
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.
155              */
156             if ((comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) && (j%2 == 0))
157             {
158                 scaling = 1.0;
159             }
160             else
161             {
162                 scaling = fudge;
163             }
164
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];
167         }
168     }
169 }
170
171 double check_mol(gmx_mtop_t *mtop, warninp_t wi)
172 {
173     char     buf[256];
174     int      i, mb, nmol, ri, pt;
175     double   q;
176     real     m;
177     t_atoms *atoms;
178
179     /* Check mass and charge */
180     q = 0.0;
181
182     for (mb = 0; mb < mtop->nmoltype; mb++)
183     {
184         atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
185         nmol  = mtop->molblock[mb].nmol;
186         for (i = 0; (i < atoms->nr); i++)
187         {
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
193              */
194             if ((m <= 0.0) && ((pt == eptAtom) || (pt == eptNucleus)))
195             {
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,
201                         m);
202                 warning_error(wi, buf);
203             }
204             else
205             if ((m != 0) && (pt == eptVSite))
206             {
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,
213                         m);
214                 warning_error(wi, buf);
215                 /* The following statements make LINCS break! */
216                 /* atoms->atom[i].m=0; */
217             }
218         }
219     }
220     return q;
221 }
222
223 static void sum_q(t_atoms *atoms, int n, double *qt, double *qBt)
224 {
225     double  qmolA, qmolB;
226     int     i;
227
228     /* sum charge */
229     qmolA = 0;
230     qmolB = 0;
231     for (i = 0; i < atoms->nr; i++)
232     {
233         qmolA += atoms->atom[i].q;
234         qmolB += atoms->atom[i].qB;
235     }
236     /* Unfortunately an absolute comparison,
237      * but this avoids unnecessary warnings and gmx-users mails.
238      */
239     if (fabs(qmolA) >= 1e-6 || fabs(qmolB) >= 1e-6)
240     {
241         *qt  += n*qmolA;
242         *qBt += n*qmolB;
243     }
244 }
245
246 static void get_nbparm(char *nb_str, char *comb_str, int *nb, int *comb,
247                        warninp_t wi)
248 {
249     int  i;
250     char warn_buf[STRLEN];
251
252     *nb   = -1;
253     for (i = 1; (i < eNBF_NR); i++)
254     {
255         if (gmx_strcasecmp(nb_str, enbf_names[i]) == 0)
256         {
257             *nb = i;
258         }
259     }
260     if (*nb == -1)
261     {
262         *nb = strtol(nb_str, NULL, 10);
263     }
264     if ((*nb < 1) || (*nb >= eNBF_NR))
265     {
266         sprintf(warn_buf, "Invalid nonbond function selector '%s' using %s",
267                 nb_str, enbf_names[1]);
268         warning_error(wi, warn_buf);
269         *nb = 1;
270     }
271     *comb = -1;
272     for (i = 1; (i < eCOMB_NR); i++)
273     {
274         if (gmx_strcasecmp(comb_str, ecomb_names[i]) == 0)
275         {
276             *comb = i;
277         }
278     }
279     if (*comb == -1)
280     {
281         *comb = strtol(comb_str, NULL, 10);
282     }
283     if ((*comb < 1) || (*comb >= eCOMB_NR))
284     {
285         sprintf(warn_buf, "Invalid combination rule selector '%s' using %s",
286                 comb_str, ecomb_names[1]);
287         warning_error(wi, warn_buf);
288         *comb = 1;
289     }
290 }
291
292 static char ** cpp_opts(const char *define, const char *include,
293                         warninp_t wi)
294 {
295     int         n, len;
296     int         ncppopts = 0;
297     const char *cppadds[2];
298     char      **cppopts   = NULL;
299     const char *option[2] = { "-D", "-I" };
300     const char *nopt[2]   = { "define", "include" };
301     const char *ptr;
302     const char *rptr;
303     char       *buf;
304     char        warn_buf[STRLEN];
305
306     cppadds[0] = define;
307     cppadds[1] = include;
308     for (n = 0; (n < 2); n++)
309     {
310         if (cppadds[n])
311         {
312             ptr = cppadds[n];
313             while (*ptr != '\0')
314             {
315                 while ((*ptr != '\0') && isspace(*ptr))
316                 {
317                     ptr++;
318                 }
319                 rptr = ptr;
320                 while ((*rptr != '\0') && !isspace(*rptr))
321                 {
322                     rptr++;
323                 }
324                 len = (rptr - ptr);
325                 if (len > 2)
326                 {
327                     snew(buf, (len+1));
328                     strncpy(buf, ptr, len);
329                     if (strstr(ptr, option[n]) != ptr)
330                     {
331                         set_warning_line(wi, "mdp file", -1);
332                         sprintf(warn_buf, "Malformed %s option %s", nopt[n], buf);
333                         warning(wi, warn_buf);
334                     }
335                     else
336                     {
337                         srenew(cppopts, ++ncppopts);
338                         cppopts[ncppopts-1] = strdup(buf);
339                     }
340                     sfree(buf);
341                     ptr = rptr;
342                 }
343             }
344         }
345     }
346     srenew(cppopts, ++ncppopts);
347     cppopts[ncppopts-1] = NULL;
348
349     return cppopts;
350 }
351
352
353 int
354 find_gb_bondlength(t_params *plist, int ai, int aj, real *length)
355 {
356     int i, j, a1, a2;
357
358     int found = 0;
359     int status;
360
361     for (i = 0; i < F_NRE && !found; i++)
362     {
363         if (IS_CHEMBOND(i))
364         {
365             for (j = 0; j < plist[i].nr; j++)
366             {
367                 a1 = plist[i].param[j].a[0];
368                 a2 = plist[i].param[j].a[1];
369
370                 if ( (a1 == ai && a2 == aj) || (a1 == aj && a2 == ai))
371                 {
372                     /* Equilibrium bond distance */
373                     *length = plist[i].param[j].c[0];
374                     found   = 1;
375                 }
376             }
377         }
378     }
379     status = !found;
380
381     return status;
382 }
383
384
385 int
386 find_gb_anglelength(t_params *plist, int ai, int ak, real *length)
387 {
388     int  i, j, a1, a2, a3;
389     real r12, r23, a123;
390     int  found = 0;
391     int  status, status1, status2;
392
393     r12 = r23 = 0;
394
395     for (i = 0; i < F_NRE && !found; i++)
396     {
397         if (IS_ANGLE(i))
398         {
399             for (j = 0; j < plist[i].nr; j++)
400             {
401                 a1 = plist[i].param[j].a[0];
402                 a2 = plist[i].param[j].a[1];
403                 a3 = plist[i].param[j].a[2];
404
405                 /* We dont care what the middle atom is, but use it below */
406                 if ( (a1 == ai && a3 == ak) || (a1 == ak && a3 == ai) )
407                 {
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);
413
414                     if (status1 == 0 && status2 == 0)
415                     {
416                         /* cosine theorem to get r13 */
417                         *length = sqrt(r12*r12+r23*r23-(2*r12*r23*cos(a123/RAD2DEG)));
418                         found   = 1;
419                     }
420                 }
421             }
422         }
423     }
424     status = !found;
425
426     return status;
427 }
428
429 int
430 generate_gb_exclusion_interactions(t_molinfo *mi, gpp_atomtype_t atype, t_nextnb *nnb)
431 {
432     int          i, j, k, n, ai, aj, ti, tj;
433     int          n12, n13, n14;
434     int          ftype;
435     t_param      param;
436     t_params *   plist;
437     t_atoms *    at;
438     real         radiusi, radiusj;
439     real         gb_radiusi, gb_radiusj;
440     real         param_c2, param_c4;
441     real         distance;
442
443     plist = mi->plist;
444     at    = &mi->atoms;
445
446     for (n = 1; n <= nnb->nrex; n++)
447     {
448         switch (n)
449         {
450             case 1:
451                 ftype    = F_GB12;
452                 param_c2 = STILL_P2;
453                 param_c4 = 0.8875;
454                 break;
455             case 2:
456                 ftype    = F_GB13;
457                 param_c2 = STILL_P3;
458                 param_c4 = 0.3516;
459                 break;
460             default:
461                 /* Put all higher-order exclusions into 1,4 list so we dont miss them */
462                 ftype    = F_GB14;
463                 param_c2 = STILL_P3;
464                 param_c4 = 0.3516;
465                 break;
466         }
467
468         for (ai = 0; ai < nnb->nr; ai++)
469         {
470             ti         = at->atom[ai].type;
471             radiusi    = get_atomtype_radius(ti, atype);
472             gb_radiusi = get_atomtype_gb_radius(ti, atype);
473
474             for (j = 0; j < nnb->nrexcl[ai][n]; j++)
475             {
476                 aj = nnb->a[ai][n][j];
477
478                 /* Only add the interactions once */
479                 if (aj > ai)
480                 {
481                     tj         = at->atom[aj].type;
482                     radiusj    = get_atomtype_radius(tj, atype);
483                     gb_radiusj = get_atomtype_gb_radius(tj, atype);
484
485                     /* There is an exclusion of type "ftype" between atoms ai and aj */
486                     param.a[0] = ai;
487                     param.a[1] = aj;
488
489                     /* Reference distance, not used for 1-4 interactions */
490                     switch (ftype)
491                     {
492                         case F_GB12:
493                             if (find_gb_bondlength(plist, ai, aj, &distance) != 0)
494                             {
495                                 gmx_fatal(FARGS, "Cannot find bond length for atoms %d-%d", ai, aj);
496                             }
497                             break;
498                         case F_GB13:
499                             if (find_gb_anglelength(plist, ai, aj, &distance) != 0)
500                             {
501                                 gmx_fatal(FARGS, "Cannot find length for atoms %d-%d involved in angle", ai, aj);
502                             }
503                             break;
504                         default:
505                             distance = -1;
506                             break;
507                     }
508                     /* Assign GB parameters */
509                     /* Sum of radii */
510                     param.c[0] = radiusi+radiusj;
511                     /* Reference distance distance */
512                     param.c[1] = distance;
513                     /* Still parameter */
514                     param.c[2] = param_c2;
515                     /* GB radius */
516                     param.c[3] = gb_radiusi+gb_radiusj;
517                     /* Parameter */
518                     param.c[4] = param_c4;
519
520                     /* Add it to the parameter list */
521                     add_param_to_list(&plist[ftype], &param);
522                 }
523             }
524         }
525     }
526     return 0;
527 }
528
529
530
531
532 static char **read_topol(const char *infile, const char *outfile,
533                          const char *define, const char *include,
534                          t_symtab    *symtab,
535                          gpp_atomtype_t atype,
536                          int         *nrmols,
537                          t_molinfo   **molinfo,
538                          t_params    plist[],
539                          int         *combination_rule,
540                          double      *reppow,
541                          t_gromppopts *opts,
542                          real        *fudgeQQ,
543                          int         *nmolblock,
544                          gmx_molblock_t **molblock,
545                          gmx_bool        bFEP,
546                          gmx_bool        bGenborn,
547                          gmx_bool        bZero,
548                          warninp_t   wi)
549 {
550     FILE           *out;
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];
554     char            genpairs[32];
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;
561     DirStack       *DS;
562     directive       d, newd;
563     t_nbparam     **nbparam, **pair;
564     t_block2       *block2;
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;
569     int             lastcg = -1;
570     int             dcatt  = -1, nmol_couple;
571     /* File handling variables */
572     int             status, done;
573     gmx_cpp_t       handle;
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.
582      */
583     if (outfile)
584     {
585         out = gmx_fio_fopen(outfile, "w");
586     }
587     else
588     {
589         out = NULL;
590     }
591
592     /* open input file */
593     status = cpp_open_file(infile, &handle, cpp_opts(define, include, wi));
594     if (status != 0)
595     {
596         gmx_fatal(FARGS, cpp_error(&handle, status));
597     }
598
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                        */
606     nb_funct = F_LJ;
607     *reppow  = 12.0;      /* Default value for repulsion power     */
608
609     comb     = 0;
610
611     /* Init the number of CMAP torsion angles  and grid spacing */
612     plist->grid_spacing = 0;
613     plist->nc           = 0;
614
615     bWarn_copy_A_B = bFEP;
616
617     batype = init_bond_atomtype();
618     /* parse the actual file */
619     bReadDefaults = FALSE;
620     bGenPairs     = FALSE;
621     bReadMolType  = FALSE;
622     nmol_couple   = 0;
623
624     do
625     {
626         status = cpp_read_line(&handle, STRLEN, line);
627         done   = (status == eCPP_EOF);
628         if (!done)
629         {
630             if (status != eCPP_OK)
631             {
632                 gmx_fatal(FARGS, cpp_error(&handle, status));
633             }
634             else if (out)
635             {
636                 fprintf(out, "%s\n", line);
637             }
638
639             set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
640
641             pline = strdup(line);
642
643             /* Strip trailing '\' from pline, if it exists */
644             sl = strlen(pline);
645             if ((sl > 0) && (pline[sl-1] == CONTINUE))
646             {
647                 pline[sl-1] = ' ';
648             }
649
650             /* build one long line from several fragments - necessary for CMAP */
651             while (continuing(line))
652             {
653                 status = cpp_read_line(&handle, STRLEN, line);
654                 set_warning_line(wi, cpp_cur_file(&handle), cpp_cur_linenr(&handle));
655
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
658                  */
659                 tmp_line = strdup(line);
660
661                 sl = strlen(tmp_line);
662                 if ((sl > 0) && (tmp_line[sl-1] == CONTINUE))
663                 {
664                     tmp_line[sl-1] = ' ';
665                 }
666
667                 done = (status == eCPP_EOF);
668                 if (!done)
669                 {
670                     if (status != eCPP_OK)
671                     {
672                         gmx_fatal(FARGS, cpp_error(&handle, status));
673                     }
674                     else if (out)
675                     {
676                         fprintf(out, "%s\n", line);
677                     }
678                 }
679
680                 srenew(pline, strlen(pline)+strlen(tmp_line)+1);
681                 strcat(pline, tmp_line);
682                 sfree(tmp_line);
683             }
684
685             /* skip trailing and leading spaces and comment text */
686             strip_comment (pline);
687             trim (pline);
688
689             /* if there is something left... */
690             if ((int)strlen(pline) > 0)
691             {
692                 if (pline[0] == OPENDIR)
693                 {
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
697                      */
698                     dirstr = strdup((pline+1));
699                     if ((dummy2 = strchr (dirstr, CLOSEDIR)) != NULL)
700                     {
701                         (*dummy2) = 0;
702                     }
703                     trim (dirstr);
704
705                     if ((newd = str2dir(dirstr)) == d_invalid)
706                     {
707                         sprintf(errbuf, "Invalid directive %s", dirstr);
708                         warning_error(wi, errbuf);
709                     }
710                     else
711                     {
712                         /* Directive found */
713                         if (debug)
714                         {
715                             fprintf(debug, "found directive '%s'\n", dir2str(newd));
716                         }
717                         if (DS_Check_Order (DS, newd))
718                         {
719                             DS_Push (&DS, newd);
720                             d = newd;
721                         }
722                         else
723                         {
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));
728                             /* d = d_invalid; */
729                         }
730                     }
731                     sfree(dirstr);
732                 }
733                 else if (d != d_invalid)
734                 {
735                     /* Not a directive, just a plain string
736                      * use a gigantic switch to decode,
737                      * if there is a valid directive!
738                      */
739                     switch (d)
740                     {
741                         case d_defaults:
742                             if (bReadDefaults)
743                             {
744                                 gmx_fatal(FARGS, "%s\nFound a second defaults directive.\n",
745                                           cpp_error(&handle, eCPP_SYNTAX));
746                             }
747                             bReadDefaults = TRUE;
748                             nscan         = sscanf(pline, "%s%s%s%lf%lf%lf",
749                                                    nb_str, comb_str, genpairs, &fLJ, &fQQ, &fPOW);
750                             if (nscan < 2)
751                             {
752                                 too_few(wi);
753                             }
754                             else
755                             {
756                                 bGenPairs = FALSE;
757                                 fudgeLJ   = 1.0;
758                                 *fudgeQQ  = 1.0;
759
760                                 get_nbparm(nb_str, comb_str, &nb_funct, &comb, wi);
761                                 *combination_rule = comb;
762                                 if (nscan >= 3)
763                                 {
764                                     bGenPairs = (gmx_strncasecmp(genpairs, "Y", 1) == 0);
765                                     if (nb_funct != eNBF_LJ && bGenPairs)
766                                     {
767                                         gmx_fatal(FARGS, "Generating pair parameters is only supported with LJ non-bonded interactions");
768                                     }
769                                 }
770                                 if (nscan >= 4)
771                                 {
772                                     fudgeLJ   = fLJ;
773                                 }
774                                 if (nscan >= 5)
775                                 {
776                                     *fudgeQQ  = fQQ;
777                                 }
778                                 if (nscan >= 6)
779                                 {
780                                     *reppow   = fPOW;
781                                 }
782                             }
783                             nb_funct = ifunc_index(d_nonbond_params, nb_funct);
784
785                             break;
786                         case d_atomtypes:
787                             push_at(symtab, atype, batype, pline, nb_funct,
788                                     &nbparam, bGenPairs ? &pair : NULL, wi);
789                             break;
790
791                         case d_bondtypes:
792                             push_bt(d, plist, 2, NULL, batype, pline, wi);
793                             break;
794                         case d_constrainttypes:
795                             push_bt(d, plist, 2, NULL, batype, pline, wi);
796                             break;
797                         case d_pairtypes:
798                             if (bGenPairs)
799                             {
800                                 push_nbt(d, pair, atype, pline, F_LJ14, wi);
801                             }
802                             else
803                             {
804                                 push_bt(d, plist, 2, atype, NULL, pline, wi);
805                             }
806                             break;
807                         case d_angletypes:
808                             push_bt(d, plist, 3, NULL, batype, pline, wi);
809                             break;
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);
813                             break;
814
815                         case d_nonbond_params:
816                             push_nbt(d, nbparam, atype, pline, nb_funct, wi);
817                             break;
818                         /*
819                            case d_blocktype:
820                            nblock++;
821                            srenew(block,nblock);
822                            srenew(blockinfo,nblock);
823                            blk0=&(block[nblock-1]);
824                            bi0=&(blockinfo[nblock-1]);
825                            init_top(blk0);
826                            init_molinfo(bi0);
827                            push_molt(symtab,bi0,pline);
828                            break;
829                          */
830
831                         case d_implicit_genborn_params:
832                             push_gb_params(atype, pline, wi);
833                             break;
834
835                         case d_implicit_surface_params:
836                             gmx_fatal(FARGS, "Implicit surface directive not supported yet.");
837                             break;
838
839                         case d_cmaptypes:
840                             push_cmaptype(d, plist, 5, atype, batype, pline, wi);
841                             break;
842
843                         case d_moleculetype:
844                         {
845                             if (!bReadMolType)
846                             {
847                                 int ntype;
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))
853                                 {
854                                     dcatt = add_atomtype_decoupled(symtab, atype,
855                                                                    &nbparam, bGenPairs ? &pair : NULL);
856                                 }
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]),
861                                                       ntype);
862                                 fprintf(stderr, "Generated %d of the %d non-bonded parameter combinations\n", ncombs-ncopy, ncombs);
863                                 free_nbparam(nbparam, ntype);
864                                 if (bGenPairs)
865                                 {
866                                     gen_pairs(&(plist[nb_funct]), &(plist[F_LJ14]), fudgeLJ, comb);
867                                     ncopy = copy_nbparams(pair, nb_funct, &(plist[F_LJ14]),
868                                                           ntype);
869                                     fprintf(stderr, "Generated %d of the %d 1-4 parameter combinations\n", ncombs-ncopy, ncombs);
870                                     free_nbparam(pair, ntype);
871                                 }
872                                 /* Copy GBSA parameters to atomtype array? */
873
874                                 bReadMolType = TRUE;
875                             }
876
877                             push_molt(symtab, &nmol, molinfo, pline, wi);
878                             srenew(block2, nmol);
879                             block2[nmol-1].nr = 0;
880                             mi0               = &((*molinfo)[nmol-1]);
881                             break;
882                         }
883                         case d_atoms:
884                             push_atom(symtab, &(mi0->cgs), &(mi0->atoms), atype, pline, &lastcg, wi);
885                             break;
886
887                         case d_pairs:
888                             push_bond(d, plist, mi0->plist, &(mi0->atoms), atype, pline, FALSE,
889                                       bGenPairs, *fudgeQQ, bZero, &bWarn_copy_A_B, wi);
890                             break;
891
892                         case d_vsites2:
893                         case d_vsites3:
894                         case d_vsites4:
895                         case d_bonds:
896                         case d_angles:
897                         case d_constraints:
898                         case d_settles:
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:
905                         case d_dihedrals:
906                         case d_polarization:
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);
911                             break;
912                         case d_cmap:
913                             push_cmap(d, plist, mi0->plist, &(mi0->atoms), atype, pline, wi);
914                             break;
915
916                         case d_vsitesn:
917                             push_vsitesn(d, mi0->plist, &(mi0->atoms), pline, wi);
918                             break;
919                         case d_exclusions:
920                             assert(block2);
921                             if (!block2[nmol-1].nr)
922                             {
923                                 init_block2(&(block2[nmol-1]), mi0->atoms.nr);
924                             }
925                             push_excl(pline, &(block2[nmol-1]));
926                             break;
927                         case d_system:
928                             trim(pline);
929                             title = put_symtab(symtab, pline);
930                             break;
931                         case d_molecules:
932                         {
933                             int      whichmol;
934                             gmx_bool bCouple;
935
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;
941                             nmolb++;
942
943                             bCouple = (opts->couple_moltype != NULL &&
944                                        (gmx_strcasecmp("system", opts->couple_moltype) == 0 ||
945                                         gmx_strcasecmp(*(mi0->name), opts->couple_moltype) == 0));
946                             if (bCouple)
947                             {
948                                 nmol_couple += nrcopies;
949                             }
950
951                             if (mi0->atoms.nr == 0)
952                             {
953                                 gmx_fatal(FARGS, "Molecule type '%s' contains no atoms",
954                                           *mi0->name);
955                             }
956                             fprintf(stderr,
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)
961                             {
962                                 t_nextnb nnb;
963                                 generate_excl(mi0->nrexcl,
964                                               mi0->atoms.nr,
965                                               mi0->plist,
966                                               &nnb,
967                                               &(mi0->excls));
968                                 merge_excl(&(mi0->excls), &(block2[whichmol]));
969                                 done_block2(&(block2[whichmol]));
970                                 make_shake(mi0->plist, &mi0->atoms, opts->nshake);
971
972
973
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.
976                                  */
977                                 if (bGenborn == TRUE)
978                                 {
979                                     generate_gb_exclusion_interactions(mi0, atype, &nnb);
980                                 }
981
982                                 done_nnb(&nnb);
983
984                                 if (bCouple)
985                                 {
986                                     convert_moltype_couple(mi0, dcatt, *fudgeQQ,
987                                                            opts->couple_lam0, opts->couple_lam1,
988                                                            opts->bCoupleIntra,
989                                                            nb_funct, &(plist[nb_funct]));
990                                 }
991                                 stupid_fill_block(&mi0->mols, mi0->atoms.nr, TRUE);
992                                 mi0->bProcessed = TRUE;
993                             }
994                             break;
995                         }
996                         default:
997                             fprintf (stderr, "case: %d\n", d);
998                             gmx_incons("unknown directive");
999                     }
1000                 }
1001             }
1002             sfree(pline);
1003             pline = NULL;
1004         }
1005     }
1006     while (!done);
1007     status = cpp_close_file(&handle);
1008     if (status != eCPP_OK)
1009     {
1010         gmx_fatal(FARGS, cpp_error(&handle, status));
1011     }
1012     cpp_done();
1013     if (out)
1014     {
1015         gmx_fio_fclose(out);
1016     }
1017
1018     if (opts->couple_moltype)
1019     {
1020         if (nmol_couple == 0)
1021         {
1022             gmx_fatal(FARGS, "Did not find any molecules of type '%s' for coupling",
1023                       opts->couple_moltype);
1024         }
1025         fprintf(stderr, "Coupling %d copies of molecule type '%s'\n",
1026                 nmol_couple, opts->couple_moltype);
1027     }
1028
1029     /* this is not very clean, but fixes core dump on empty system name */
1030     if (!title)
1031     {
1032         title = put_symtab(symtab, "");
1033     }
1034     if (fabs(qt) > 1e-4)
1035     {
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);
1038     }
1039     if (fabs(qBt) > 1e-4 && !gmx_within_tol(qBt, qt, 1e-6))
1040     {
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);
1043     }
1044     DS_Done (&DS);
1045     for (i = 0; i < nmol; i++)
1046     {
1047         done_block2(&(block2[i]));
1048     }
1049     free(block2);
1050
1051     done_bond_atomtype(&batype);
1052
1053     *nrmols = nmol;
1054
1055     *nmolblock = nmolb;
1056     *molblock  = molb;
1057
1058     return title;
1059 }
1060
1061 char **do_top(gmx_bool          bVerbose,
1062               const char       *topfile,
1063               const char       *topppfile,
1064               t_gromppopts     *opts,
1065               gmx_bool          bZero,
1066               t_symtab         *symtab,
1067               t_params          plist[],
1068               int              *combination_rule,
1069               double           *repulsion_power,
1070               real             *fudgeQQ,
1071               gpp_atomtype_t    atype,
1072               int              *nrmols,
1073               t_molinfo       **molinfo,
1074               t_inputrec       *ir,
1075               int              *nmolblock,
1076               gmx_molblock_t  **molblock,
1077               gmx_bool          bGenborn,
1078               warninp_t         wi)
1079 {
1080     /* Tmpfile might contain a long path */
1081     const char *tmpfile;
1082     char      **title;
1083
1084     if (topppfile)
1085     {
1086         tmpfile = topppfile;
1087     }
1088     else
1089     {
1090         tmpfile = NULL;
1091     }
1092
1093     if (bVerbose)
1094     {
1095         printf("processing topology...\n");
1096     }
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))
1104     {
1105         warning(wi, "Using sigma/epsilon based combination rules with"
1106                 " user supplied potential function may produce unwanted"
1107                 " results");
1108     }
1109
1110     return title;
1111 }
1112
1113
1114 static void generate_qmexcl_moltype(gmx_moltype_t *molt, unsigned char *grpnr,
1115                                     t_inputrec *ir)
1116 {
1117     /* This routine expects molt->ilist to be of size F_NRE and ordered. */
1118
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
1122      */
1123     int
1124         i, j, l, k = 0, jmax, qm_max = 0, qm_nr = 0, nratoms = 0, link_nr = 0, link_max = 0;
1125     atom_id
1126        *qm_arr = NULL, *link_arr = NULL, a1, a2, a3, a4, ftype = 0;
1127     t_blocka
1128         qmexcl;
1129     t_block2
1130         qmexcl2;
1131     gmx_bool
1132        *bQMMM, *blink, bexcl;
1133
1134     /* First we search and select the QM atoms in an qm_arr array that
1135      * we use to create the exclusions.
1136      *
1137      * we take the possibility into account that a user has defined more
1138      * than one QM group:
1139      *
1140      * for that we also need to do this an ugly work-about just in case
1141      * the QM group contains the entire system...
1142      */
1143     jmax = ir->opts.ngQM;
1144
1145     /* we first search for all the QM atoms and put them in an array
1146      */
1147     for (j = 0; j < jmax; j++)
1148     {
1149         for (i = 0; i < molt->atoms.nr; i++)
1150         {
1151             if (qm_nr >= qm_max)
1152             {
1153                 qm_max += 100;
1154                 srenew(qm_arr, qm_max);
1155             }
1156             if ((grpnr ? grpnr[i] : 0) == j)
1157             {
1158                 qm_arr[qm_nr++] = i;
1159             }
1160         }
1161     }
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
1165      * to TRUE.
1166      */
1167     snew(bQMMM, molt->atoms.nr);
1168     for (i = 0; i < molt->atoms.nr; i++)
1169     {
1170         bQMMM[i] = FALSE;
1171     }
1172     for (i = 0; i < qm_nr; i++)
1173     {
1174         bQMMM[qm_arr[i]] = TRUE;
1175     }
1176
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 */
1186
1187     /* first check weter we already have CONNBONDS: */
1188     if (molt->ilist[F_CONNBONDS].nr != 0)
1189     {
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;
1194     }
1195     /* now we delete all bonded interactions, except the ones describing
1196      * a chemical bond. These are converted to CONNBONDS
1197      */
1198     for (i = 0; i < F_LJ; i++)
1199     {
1200         if (i == F_CONNBONDS)
1201         {
1202             continue;
1203         }
1204         nratoms = interaction_function[i].nratoms;
1205         j       = 0;
1206         while (j < molt->ilist[i].nr)
1207         {
1208             bexcl = FALSE;
1209             switch (nratoms)
1210             {
1211                 case 2:
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
1217                      */
1218                     if (bexcl && i < F_ANGLES)
1219                     {
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;
1225                     }
1226                     break;
1227                 case 3:
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]));
1234                     break;
1235                 case 4:
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]));
1244                     break;
1245                 default:
1246                     gmx_fatal(FARGS, "no such bonded interactions with %d atoms\n", nratoms);
1247             }
1248             if (bexcl)
1249             {
1250                 /* since the interaction involves QM atoms, these should be
1251                  * removed from the MM ilist
1252                  */
1253                 molt->ilist[i].nr -= (nratoms+1);
1254                 for (l = j; l < molt->ilist[i].nr; l++)
1255                 {
1256                     molt->ilist[i].iatoms[l] = molt->ilist[i].iatoms[l+(nratoms+1)];
1257                 }
1258             }
1259             else
1260             {
1261                 j += nratoms+1; /* the +1 is for the functype */
1262             }
1263         }
1264     }
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
1269      * twice.  */
1270
1271     for (i = 0; i < F_NRE; i++)
1272     {
1273         if (IS_CHEMBOND(i))
1274         {
1275             j = 0;
1276             while (j < molt->ilist[i].nr)
1277             {
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]))
1281                 {
1282                     if (link_nr >= link_max)
1283                     {
1284                         link_max += 10;
1285                         srenew(link_arr, link_max);
1286                     }
1287                     if (bQMMM[a1])
1288                     {
1289                         link_arr[link_nr++] = a2;
1290                     }
1291                     else
1292                     {
1293                         link_arr[link_nr++] = a1;
1294                     }
1295                 }
1296                 j += 3;
1297             }
1298         }
1299     }
1300     snew(blink, molt->atoms.nr);
1301     for (i = 0; i < molt->atoms.nr; i++)
1302     {
1303         blink[i] = FALSE;
1304     }
1305     for (i = 0; i < link_nr; i++)
1306     {
1307         blink[link_arr[i]] = TRUE;
1308     }
1309     /* creating the exclusion block for the QM atoms. Each QM atom has
1310      * as excluded elements all the other QMatoms (and itself).
1311      */
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);
1316     j = 0;
1317     for (i = 0; i < qmexcl.nr; i++)
1318     {
1319         qmexcl.index[i] = j;
1320         if (bQMMM[i])
1321         {
1322             for (k = 0; k < qm_nr; k++)
1323             {
1324                 qmexcl.a[k+j] = qm_arr[k];
1325             }
1326             for (k = 0; k < link_nr; k++)
1327             {
1328                 qmexcl.a[qm_nr+k+j] = link_arr[k];
1329             }
1330             j += (qm_nr+link_nr);
1331         }
1332         if (blink[i])
1333         {
1334             for (k = 0; k < qm_nr; k++)
1335             {
1336                 qmexcl.a[k+j] = qm_arr[k];
1337             }
1338             j += qm_nr;
1339         }
1340     }
1341     qmexcl.index[qmexcl.nr] = j;
1342
1343     /* and merging with the exclusions already present in sys.
1344      */
1345
1346     init_block2(&qmexcl2, molt->atoms.nr);
1347     b_to_b2(&qmexcl, &qmexcl2);
1348     merge_excl(&(molt->excls), &qmexcl2);
1349     done_block2(&qmexcl2);
1350
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++)
1357     {
1358         nratoms = interaction_function[i].nratoms;
1359         j       = 0;
1360         while (j < molt->ilist[i].nr)
1361         {
1362             bexcl = FALSE;
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]));
1368             if (bexcl)
1369             {
1370                 /* since the interaction involves QM atoms, these should be
1371                  * removed from the MM ilist
1372                  */
1373                 molt->ilist[i].nr -= (nratoms+1);
1374                 for (k = j; k < molt->ilist[i].nr; k++)
1375                 {
1376                     molt->ilist[i].iatoms[k] = molt->ilist[i].iatoms[k+(nratoms+1)];
1377                 }
1378             }
1379             else
1380             {
1381                 j += nratoms+1; /* the +1 is for the functype */
1382             }
1383         }
1384     }
1385
1386     free(qm_arr);
1387     free(bQMMM);
1388     free(link_arr);
1389     free(blink);
1390 } /* generate_qmexcl */
1391
1392 void generate_qmexcl(gmx_mtop_t *sys, t_inputrec *ir, warninp_t    wi)
1393 {
1394     /* This routine expects molt->molt[m].ilist to be of size F_NRE and ordered.
1395      */
1396
1397     unsigned char  *grpnr;
1398     int             mb, mol, nat_mol, i, nr_mol_with_qm_atoms = 0;
1399     gmx_molblock_t *molb;
1400     gmx_bool        bQMMM;
1401
1402     grpnr = sys->groups.grpnr[egcQMMM];
1403
1404     for (mb = 0; mb < sys->nmolblock; mb++)
1405     {
1406         molb    = &sys->molblock[mb];
1407         nat_mol = sys->moltype[molb->type].atoms.nr;
1408         for (mol = 0; mol < molb->nmol; mol++)
1409         {
1410             bQMMM = FALSE;
1411             for (i = 0; i < nat_mol; i++)
1412             {
1413                 if ((grpnr ? grpnr[i] : 0) < ir->opts.ngQM)
1414                 {
1415                     bQMMM = TRUE;
1416                 }
1417             }
1418             if (bQMMM)
1419             {
1420                 nr_mol_with_qm_atoms++;
1421                 if (molb->nmol > 1)
1422                 {
1423                     /* We need to split this molblock */
1424                     if (mol > 0)
1425                     {
1426                         /* Split the molblock at this molecule */
1427                         sys->nmolblock++;
1428                         srenew(sys->molblock, sys->nmolblock);
1429                         for (i = sys->nmolblock-2; i >= mb; i--)
1430                         {
1431                             sys->molblock[i+1] = sys->molblock[i];
1432                         }
1433                         sys->molblock[mb  ].nmol  = mol;
1434                         sys->molblock[mb+1].nmol -= mol;
1435                         mb++;
1436                         molb = &sys->molblock[mb];
1437                     }
1438                     if (molb->nmol > 1)
1439                     {
1440                         /* Split the molblock after this molecule */
1441                         sys->nmolblock++;
1442                         srenew(sys->molblock, sys->nmolblock);
1443                         molb = &sys->molblock[mb];
1444                         for (i = sys->nmolblock-2; i >= mb; i--)
1445                         {
1446                             sys->molblock[i+1] = sys->molblock[i];
1447                         }
1448                         sys->molblock[mb  ].nmol  = 1;
1449                         sys->molblock[mb+1].nmol -= 1;
1450                     }
1451
1452                     /* Add a moltype for the QMMM molecule */
1453                     sys->nmoltype++;
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.
1459                      */
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;
1464                 }
1465                 generate_qmexcl_moltype(&sys->moltype[molb->type], grpnr, ir);
1466             }
1467             if (grpnr)
1468             {
1469                 grpnr += nat_mol;
1470             }
1471         }
1472     }
1473     if (nr_mol_with_qm_atoms > 1)
1474     {
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
1480          */
1481         warning_note(wi,
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"
1491                      "\n\n");
1492     }
1493 }