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