Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / constr.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 <stdlib.h>
42
43 #include "gromacs/fileio/confio.h"
44 #include "types/commrec.h"
45 #include "constr.h"
46 #include "copyrite.h"
47 #include "mdrun.h"
48 #include "nrnb.h"
49 #include "gromacs/math/vec.h"
50 #include "names.h"
51 #include "txtdump.h"
52 #include "domdec.h"
53 #include "gromacs/fileio/pdbio.h"
54 #include "splitter.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/fileio/gmxfio.h"
57 #include "macros.h"
58 #include "gmx_omp_nthreads.h"
59 #include "gromacs/essentialdynamics/edsam.h"
60 #include "gromacs/pulling/pull.h"
61
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/topology/block.h"
64 #include "gromacs/topology/invblock.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/smalloc.h"
67
68 typedef struct gmx_constr {
69     int                ncon_tot;       /* The total number of constraints    */
70     int                nflexcon;       /* The number of flexible constraints */
71     int                n_at2con_mt;    /* The size of at2con = #moltypes     */
72     t_blocka          *at2con_mt;      /* A list of atoms to constraints     */
73     int                n_at2settle_mt; /* The size of at2settle = #moltypes  */
74     int              **at2settle_mt;   /* A list of atoms to settles         */
75     gmx_bool           bInterCGsettles;
76     gmx_lincsdata_t    lincsd;         /* LINCS data                         */
77     gmx_shakedata_t    shaked;         /* SHAKE data                         */
78     gmx_settledata_t   settled;        /* SETTLE data                        */
79     int                nblocks;        /* The number of SHAKE blocks         */
80     int               *sblock;         /* The SHAKE blocks                   */
81     int                sblock_nalloc;  /* The allocation size of sblock      */
82     real              *lagr;           /* Lagrange multipliers for SHAKE     */
83     int                lagr_nalloc;    /* The allocation size of lagr        */
84     int                maxwarn;        /* The maximum number of warnings     */
85     int                warncount_lincs;
86     int                warncount_settle;
87     gmx_edsam_t        ed;            /* The essential dynamics data        */
88
89     tensor            *vir_r_m_dr_th; /* Thread local working data          */
90     int               *settle_error;  /* Thread local working data          */
91
92     gmx_mtop_t        *warn_mtop;     /* Only used for printing warnings    */
93 } t_gmx_constr;
94
95 typedef struct {
96     atom_id iatom[3];
97     atom_id blocknr;
98 } t_sortblock;
99
100 static void *init_vetavars(t_vetavars *vars,
101                            gmx_bool constr_deriv,
102                            real veta, real vetanew, t_inputrec *ir, gmx_ekindata_t *ekind, gmx_bool bPscal)
103 {
104     double g;
105     int    i;
106
107     /* first, set the alpha integrator variable */
108     if ((ir->opts.nrdf[0] > 0) && bPscal)
109     {
110         vars->alpha = 1.0 + DIM/((double)ir->opts.nrdf[0]);
111     }
112     else
113     {
114         vars->alpha = 1.0;
115     }
116     g             = 0.5*veta*ir->delta_t;
117     vars->rscale  = exp(g)*series_sinhx(g);
118     g             = -0.25*vars->alpha*veta*ir->delta_t;
119     vars->vscale  = exp(g)*series_sinhx(g);
120     vars->rvscale = vars->vscale*vars->rscale;
121     vars->veta    = vetanew;
122
123     if (constr_deriv)
124     {
125         snew(vars->vscale_nhc, ir->opts.ngtc);
126         if ((ekind == NULL) || (!bPscal))
127         {
128             for (i = 0; i < ir->opts.ngtc; i++)
129             {
130                 vars->vscale_nhc[i] = 1;
131             }
132         }
133         else
134         {
135             for (i = 0; i < ir->opts.ngtc; i++)
136             {
137                 vars->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
138             }
139         }
140     }
141     else
142     {
143         vars->vscale_nhc = NULL;
144     }
145
146     return vars;
147 }
148
149 static void free_vetavars(t_vetavars *vars)
150 {
151     if (vars->vscale_nhc != NULL)
152     {
153         sfree(vars->vscale_nhc);
154     }
155 }
156
157 static int pcomp(const void *p1, const void *p2)
158 {
159     int          db;
160     atom_id      min1, min2, max1, max2;
161     t_sortblock *a1 = (t_sortblock *)p1;
162     t_sortblock *a2 = (t_sortblock *)p2;
163
164     db = a1->blocknr-a2->blocknr;
165
166     if (db != 0)
167     {
168         return db;
169     }
170
171     min1 = min(a1->iatom[1], a1->iatom[2]);
172     max1 = max(a1->iatom[1], a1->iatom[2]);
173     min2 = min(a2->iatom[1], a2->iatom[2]);
174     max2 = max(a2->iatom[1], a2->iatom[2]);
175
176     if (min1 == min2)
177     {
178         return max1-max2;
179     }
180     else
181     {
182         return min1-min2;
183     }
184 }
185
186 int n_flexible_constraints(struct gmx_constr *constr)
187 {
188     int nflexcon;
189
190     if (constr)
191     {
192         nflexcon = constr->nflexcon;
193     }
194     else
195     {
196         nflexcon = 0;
197     }
198
199     return nflexcon;
200 }
201
202 void too_many_constraint_warnings(int eConstrAlg, int warncount)
203 {
204     const char *abort = "- aborting to avoid logfile runaway.\n"
205         "This normally happens when your system is not sufficiently equilibrated,"
206         "or if you are changing lambda too fast in free energy simulations.\n";
207
208     gmx_fatal(FARGS,
209               "Too many %s warnings (%d)\n"
210               "If you know what you are doing you can %s"
211               "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
212               "but normally it is better to fix the problem",
213               (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
214               (eConstrAlg == econtLINCS) ?
215               "adjust the lincs warning threshold in your mdp file\nor " : "\n");
216 }
217
218 static void write_constr_pdb(const char *fn, const char *title,
219                              gmx_mtop_t *mtop,
220                              int start, int homenr, t_commrec *cr,
221                              rvec x[], matrix box)
222 {
223     char          fname[STRLEN];
224     FILE         *out;
225     int           dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
226     gmx_domdec_t *dd;
227     char         *anm, *resnm;
228
229     dd = NULL;
230     if (DOMAINDECOMP(cr))
231     {
232         dd = cr->dd;
233         dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
234         start  = 0;
235         homenr = dd_ac1;
236     }
237
238     if (PAR(cr))
239     {
240         sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
241     }
242     else
243     {
244         sprintf(fname, "%s.pdb", fn);
245     }
246
247     out = gmx_fio_fopen(fname, "w");
248
249     fprintf(out, "TITLE     %s\n", title);
250     gmx_write_pdb_box(out, -1, box);
251     for (i = start; i < start+homenr; i++)
252     {
253         if (dd != NULL)
254         {
255             if (i >= dd->nat_home && i < dd_ac0)
256             {
257                 continue;
258             }
259             ii = dd->gatindex[i];
260         }
261         else
262         {
263             ii = i;
264         }
265         gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm);
266         gmx_fprintf_pdb_atomline(out, epdbATOM, ii+1, anm, ' ', resnm, ' ', resnr, ' ',
267                                  10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 0.0, "");
268     }
269     fprintf(out, "TER\n");
270
271     gmx_fio_fclose(out);
272 }
273
274 static void dump_confs(FILE *fplog, gmx_int64_t step, gmx_mtop_t *mtop,
275                        int start, int homenr, t_commrec *cr,
276                        rvec x[], rvec xprime[], matrix box)
277 {
278     char  buf[256], buf2[22];
279
280     char *env = getenv("GMX_SUPPRESS_DUMP");
281     if (env)
282     {
283         return;
284     }
285
286     sprintf(buf, "step%sb", gmx_step_str(step, buf2));
287     write_constr_pdb(buf, "initial coordinates",
288                      mtop, start, homenr, cr, x, box);
289     sprintf(buf, "step%sc", gmx_step_str(step, buf2));
290     write_constr_pdb(buf, "coordinates after constraining",
291                      mtop, start, homenr, cr, xprime, box);
292     if (fplog)
293     {
294         fprintf(fplog, "Wrote pdb files with previous and current coordinates\n");
295     }
296     fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
297 }
298
299 static void pr_sortblock(FILE *fp, const char *title, int nsb, t_sortblock sb[])
300 {
301     int i;
302
303     fprintf(fp, "%s\n", title);
304     for (i = 0; (i < nsb); i++)
305     {
306         fprintf(fp, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
307                 i, sb[i].iatom[0], sb[i].iatom[1], sb[i].iatom[2],
308                 sb[i].blocknr);
309     }
310 }
311
312 gmx_bool constrain(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
313                    struct gmx_constr *constr,
314                    t_idef *idef, t_inputrec *ir, gmx_ekindata_t *ekind,
315                    t_commrec *cr,
316                    gmx_int64_t step, int delta_step,
317                    t_mdatoms *md,
318                    rvec *x, rvec *xprime, rvec *min_proj,
319                    gmx_bool bMolPBC, matrix box,
320                    real lambda, real *dvdlambda,
321                    rvec *v, tensor *vir,
322                    t_nrnb *nrnb, int econq, gmx_bool bPscal,
323                    real veta, real vetanew)
324 {
325     gmx_bool    bOK, bDump;
326     int         start, homenr, nrend;
327     int         i, j, d;
328     int         ncons, settle_error;
329     tensor      vir_r_m_dr;
330     rvec       *vstor;
331     real        invdt, vir_fac, t;
332     t_ilist    *settle;
333     int         nsettle;
334     t_pbc       pbc, *pbc_null;
335     char        buf[22];
336     t_vetavars  vetavar;
337     int         nth, th;
338
339     if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
340     {
341         gmx_incons("constrain called for forces displacements while not doing energy minimization, can not do this while the LINCS and SETTLE constraint connection matrices are mass weighted");
342     }
343
344     bOK   = TRUE;
345     bDump = FALSE;
346
347     start  = 0;
348     homenr = md->homenr;
349     nrend  = start+homenr;
350
351     /* set constants for pressure control integration */
352     init_vetavars(&vetavar, econq != econqCoord,
353                   veta, vetanew, ir, ekind, bPscal);
354
355     if (ir->delta_t == 0)
356     {
357         invdt = 0;
358     }
359     else
360     {
361         invdt  = 1/ir->delta_t;
362     }
363
364     if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
365     {
366         /* Set the constraint lengths for the step at which this configuration
367          * is meant to be. The invmasses should not be changed.
368          */
369         lambda += delta_step*ir->fepvals->delta_lambda;
370     }
371
372     if (vir != NULL)
373     {
374         clear_mat(vir_r_m_dr);
375     }
376
377     where();
378
379     settle  = &idef->il[F_SETTLE];
380     nsettle = settle->nr/(1+NRAL(F_SETTLE));
381
382     if (nsettle > 0)
383     {
384         nth = gmx_omp_nthreads_get(emntSETTLE);
385     }
386     else
387     {
388         nth = 1;
389     }
390
391     if (nth > 1 && constr->vir_r_m_dr_th == NULL)
392     {
393         snew(constr->vir_r_m_dr_th, nth);
394         snew(constr->settle_error, nth);
395     }
396
397     settle_error = -1;
398
399     /* We do not need full pbc when constraints do not cross charge groups,
400      * i.e. when dd->constraint_comm==NULL.
401      * Note that PBC for constraints is different from PBC for bondeds.
402      * For constraints there is both forward and backward communication.
403      */
404     if (ir->ePBC != epbcNONE &&
405         (cr->dd || bMolPBC) && !(cr->dd && cr->dd->constraint_comm == NULL))
406     {
407         /* With pbc=screw the screw has been changed to a shift
408          * by the constraint coordinate communication routine,
409          * so that here we can use normal pbc.
410          */
411         pbc_null = set_pbc_dd(&pbc, ir->ePBC, cr->dd, FALSE, box);
412     }
413     else
414     {
415         pbc_null = NULL;
416     }
417
418     /* Communicate the coordinates required for the non-local constraints
419      * for LINCS and/or SETTLE.
420      */
421     if (cr->dd)
422     {
423         dd_move_x_constraints(cr->dd, box, x, xprime, econq == econqCoord);
424     }
425
426     if (constr->lincsd != NULL)
427     {
428         bOK = constrain_lincs(fplog, bLog, bEner, ir, step, constr->lincsd, md, cr,
429                               x, xprime, min_proj,
430                               box, pbc_null, lambda, dvdlambda,
431                               invdt, v, vir != NULL, vir_r_m_dr,
432                               econq, nrnb,
433                               constr->maxwarn, &constr->warncount_lincs);
434         if (!bOK && constr->maxwarn >= 0)
435         {
436             if (fplog != NULL)
437             {
438                 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
439                         econstr_names[econtLINCS], gmx_step_str(step, buf));
440             }
441             bDump = TRUE;
442         }
443     }
444
445     if (constr->nblocks > 0)
446     {
447         switch (econq)
448         {
449             case (econqCoord):
450                 bOK = bshakef(fplog, constr->shaked,
451                               md->invmass, constr->nblocks, constr->sblock,
452                               idef, ir, x, xprime, nrnb,
453                               constr->lagr, lambda, dvdlambda,
454                               invdt, v, vir != NULL, vir_r_m_dr,
455                               constr->maxwarn >= 0, econq, &vetavar);
456                 break;
457             case (econqVeloc):
458                 bOK = bshakef(fplog, constr->shaked,
459                               md->invmass, constr->nblocks, constr->sblock,
460                               idef, ir, x, min_proj, nrnb,
461                               constr->lagr, lambda, dvdlambda,
462                               invdt, NULL, vir != NULL, vir_r_m_dr,
463                               constr->maxwarn >= 0, econq, &vetavar);
464                 break;
465             default:
466                 gmx_fatal(FARGS, "Internal error, SHAKE called for constraining something else than coordinates");
467                 break;
468         }
469
470         if (!bOK && constr->maxwarn >= 0)
471         {
472             if (fplog != NULL)
473             {
474                 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
475                         econstr_names[econtSHAKE], gmx_step_str(step, buf));
476             }
477             bDump = TRUE;
478         }
479     }
480
481     if (nsettle > 0)
482     {
483         int calcvir_atom_end;
484
485         if (vir == NULL)
486         {
487             calcvir_atom_end = 0;
488         }
489         else
490         {
491             calcvir_atom_end = md->homenr;
492         }
493
494         switch (econq)
495         {
496             case econqCoord:
497 #pragma omp parallel for num_threads(nth) schedule(static)
498                 for (th = 0; th < nth; th++)
499                 {
500                     int start_th, end_th;
501
502                     if (th > 0)
503                     {
504                         clear_mat(constr->vir_r_m_dr_th[th]);
505                     }
506
507                     start_th = (nsettle* th   )/nth;
508                     end_th   = (nsettle*(th+1))/nth;
509                     if (start_th >= 0 && end_th - start_th > 0)
510                     {
511                         csettle(constr->settled,
512                                 end_th-start_th,
513                                 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
514                                 pbc_null,
515                                 x[0], xprime[0],
516                                 invdt, v ? v[0] : NULL, calcvir_atom_end,
517                                 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
518                                 th == 0 ? &settle_error : &constr->settle_error[th],
519                                 &vetavar);
520                     }
521                 }
522                 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
523                 if (v != NULL)
524                 {
525                     inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
526                 }
527                 if (vir != NULL)
528                 {
529                     inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
530                 }
531                 break;
532             case econqVeloc:
533             case econqDeriv:
534             case econqForce:
535             case econqForceDispl:
536 #pragma omp parallel for num_threads(nth) schedule(static)
537                 for (th = 0; th < nth; th++)
538                 {
539                     int start_th, end_th;
540
541                     if (th > 0)
542                     {
543                         clear_mat(constr->vir_r_m_dr_th[th]);
544                     }
545
546                     start_th = (nsettle* th   )/nth;
547                     end_th   = (nsettle*(th+1))/nth;
548
549                     if (start_th >= 0 && end_th - start_th > 0)
550                     {
551                         settle_proj(constr->settled, econq,
552                                     end_th-start_th,
553                                     settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
554                                     pbc_null,
555                                     x,
556                                     xprime, min_proj, calcvir_atom_end,
557                                     th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
558                                     &vetavar);
559                     }
560                 }
561                 /* This is an overestimate */
562                 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
563                 break;
564             case econqDeriv_FlexCon:
565                 /* Nothing to do, since the are no flexible constraints in settles */
566                 break;
567             default:
568                 gmx_incons("Unknown constraint quantity for settle");
569         }
570     }
571
572     if (settle->nr > 0)
573     {
574         /* Combine virial and error info of the other threads */
575         for (i = 1; i < nth; i++)
576         {
577             m_add(vir_r_m_dr, constr->vir_r_m_dr_th[i], vir_r_m_dr);
578             settle_error = constr->settle_error[i];
579         }
580
581         if (econq == econqCoord && settle_error >= 0)
582         {
583             bOK = FALSE;
584             if (constr->maxwarn >= 0)
585             {
586                 char buf[256];
587                 sprintf(buf,
588                         "\nstep " "%"GMX_PRId64 ": Water molecule starting at atom %d can not be "
589                         "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
590                         step, ddglatnr(cr->dd, settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
591                 if (fplog)
592                 {
593                     fprintf(fplog, "%s", buf);
594                 }
595                 fprintf(stderr, "%s", buf);
596                 constr->warncount_settle++;
597                 if (constr->warncount_settle > constr->maxwarn)
598                 {
599                     too_many_constraint_warnings(-1, constr->warncount_settle);
600                 }
601                 bDump = TRUE;
602             }
603         }
604     }
605
606     free_vetavars(&vetavar);
607
608     if (vir != NULL)
609     {
610         switch (econq)
611         {
612             case econqCoord:
613                 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
614                 break;
615             case econqVeloc:
616                 vir_fac = 0.5/ir->delta_t;
617                 break;
618             case econqForce:
619             case econqForceDispl:
620                 vir_fac = 0.5;
621                 break;
622             default:
623                 vir_fac = 0;
624                 gmx_incons("Unsupported constraint quantity for virial");
625         }
626
627         if (EI_VV(ir->eI))
628         {
629             vir_fac *= 2;  /* only constraining over half the distance here */
630         }
631         for (i = 0; i < DIM; i++)
632         {
633             for (j = 0; j < DIM; j++)
634             {
635                 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
636             }
637         }
638     }
639
640     if (bDump)
641     {
642         dump_confs(fplog, step, constr->warn_mtop, start, homenr, cr, x, xprime, box);
643     }
644
645     if (econq == econqCoord)
646     {
647         if (ir->ePull == epullCONSTRAINT)
648         {
649             if (EI_DYNAMICS(ir->eI))
650             {
651                 t = ir->init_t + (step + delta_step)*ir->delta_t;
652             }
653             else
654             {
655                 t = ir->init_t;
656             }
657             set_pbc(&pbc, ir->ePBC, box);
658             pull_constraint(ir->pull, md, &pbc, cr, ir->delta_t, t, x, xprime, v, *vir);
659         }
660         if (constr->ed && delta_step > 0)
661         {
662             /* apply the essential dynamcs constraints here */
663             do_edsam(ir, step, cr, xprime, v, box, constr->ed);
664         }
665     }
666
667     return bOK;
668 }
669
670 real *constr_rmsd_data(struct gmx_constr *constr)
671 {
672     if (constr->lincsd)
673     {
674         return lincs_rmsd_data(constr->lincsd);
675     }
676     else
677     {
678         return NULL;
679     }
680 }
681
682 real constr_rmsd(struct gmx_constr *constr, gmx_bool bSD2)
683 {
684     if (constr->lincsd)
685     {
686         return lincs_rmsd(constr->lincsd, bSD2);
687     }
688     else
689     {
690         return 0;
691     }
692 }
693
694 static void make_shake_sblock_serial(struct gmx_constr *constr,
695                                      t_idef *idef, t_mdatoms *md)
696 {
697     int          i, j, m, ncons;
698     int          bstart, bnr;
699     t_blocka     sblocks;
700     t_sortblock *sb;
701     t_iatom     *iatom;
702     atom_id     *inv_sblock;
703
704     /* Since we are processing the local topology,
705      * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
706      */
707     ncons = idef->il[F_CONSTR].nr/3;
708
709     init_blocka(&sblocks);
710     gen_sblocks(NULL, 0, md->homenr, idef, &sblocks, FALSE);
711
712     /*
713        bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
714        nblocks=blocks->multinr[idef->nodeid] - bstart;
715      */
716     bstart          = 0;
717     constr->nblocks = sblocks.nr;
718     if (debug)
719     {
720         fprintf(debug, "ncons: %d, bstart: %d, nblocks: %d\n",
721                 ncons, bstart, constr->nblocks);
722     }
723
724     /* Calculate block number for each atom */
725     inv_sblock = make_invblocka(&sblocks, md->nr);
726
727     done_blocka(&sblocks);
728
729     /* Store the block number in temp array and
730      * sort the constraints in order of the sblock number
731      * and the atom numbers, really sorting a segment of the array!
732      */
733 #ifdef DEBUGIDEF
734     pr_idef(fplog, 0, "Before Sort", idef);
735 #endif
736     iatom = idef->il[F_CONSTR].iatoms;
737     snew(sb, ncons);
738     for (i = 0; (i < ncons); i++, iatom += 3)
739     {
740         for (m = 0; (m < 3); m++)
741         {
742             sb[i].iatom[m] = iatom[m];
743         }
744         sb[i].blocknr = inv_sblock[iatom[1]];
745     }
746
747     /* Now sort the blocks */
748     if (debug)
749     {
750         pr_sortblock(debug, "Before sorting", ncons, sb);
751         fprintf(debug, "Going to sort constraints\n");
752     }
753
754     qsort(sb, ncons, (size_t)sizeof(*sb), pcomp);
755
756     if (debug)
757     {
758         pr_sortblock(debug, "After sorting", ncons, sb);
759     }
760
761     iatom = idef->il[F_CONSTR].iatoms;
762     for (i = 0; (i < ncons); i++, iatom += 3)
763     {
764         for (m = 0; (m < 3); m++)
765         {
766             iatom[m] = sb[i].iatom[m];
767         }
768     }
769 #ifdef DEBUGIDEF
770     pr_idef(fplog, 0, "After Sort", idef);
771 #endif
772
773     j = 0;
774     snew(constr->sblock, constr->nblocks+1);
775     bnr = -2;
776     for (i = 0; (i < ncons); i++)
777     {
778         if (sb[i].blocknr != bnr)
779         {
780             bnr                 = sb[i].blocknr;
781             constr->sblock[j++] = 3*i;
782         }
783     }
784     /* Last block... */
785     constr->sblock[j++] = 3*ncons;
786
787     if (j != (constr->nblocks+1))
788     {
789         fprintf(stderr, "bstart: %d\n", bstart);
790         fprintf(stderr, "j: %d, nblocks: %d, ncons: %d\n",
791                 j, constr->nblocks, ncons);
792         for (i = 0; (i < ncons); i++)
793         {
794             fprintf(stderr, "i: %5d  sb[i].blocknr: %5d\n", i, sb[i].blocknr);
795         }
796         for (j = 0; (j <= constr->nblocks); j++)
797         {
798             fprintf(stderr, "sblock[%3d]=%5d\n", j, (int)constr->sblock[j]);
799         }
800         gmx_fatal(FARGS, "DEATH HORROR: "
801                   "sblocks does not match idef->il[F_CONSTR]");
802     }
803     sfree(sb);
804     sfree(inv_sblock);
805 }
806
807 static void make_shake_sblock_dd(struct gmx_constr *constr,
808                                  t_ilist *ilcon, t_block *cgs,
809                                  gmx_domdec_t *dd)
810 {
811     int      ncons, c, cg;
812     t_iatom *iatom;
813
814     if (dd->ncg_home+1 > constr->sblock_nalloc)
815     {
816         constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
817         srenew(constr->sblock, constr->sblock_nalloc);
818     }
819
820     ncons           = ilcon->nr/3;
821     iatom           = ilcon->iatoms;
822     constr->nblocks = 0;
823     cg              = 0;
824     for (c = 0; c < ncons; c++)
825     {
826         if (c == 0 || iatom[1] >= cgs->index[cg+1])
827         {
828             constr->sblock[constr->nblocks++] = 3*c;
829             while (iatom[1] >= cgs->index[cg+1])
830             {
831                 cg++;
832             }
833         }
834         iatom += 3;
835     }
836     constr->sblock[constr->nblocks] = 3*ncons;
837 }
838
839 t_blocka make_at2con(int start, int natoms,
840                      t_ilist *ilist, t_iparams *iparams,
841                      gmx_bool bDynamics, int *nflexiblecons)
842 {
843     int      *count, ncon, con, con_tot, nflexcon, ftype, i, a;
844     t_iatom  *ia;
845     t_blocka  at2con;
846     gmx_bool  bFlexCon;
847
848     snew(count, natoms);
849     nflexcon = 0;
850     for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
851     {
852         ncon = ilist[ftype].nr/3;
853         ia   = ilist[ftype].iatoms;
854         for (con = 0; con < ncon; con++)
855         {
856             bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
857                         iparams[ia[0]].constr.dB == 0);
858             if (bFlexCon)
859             {
860                 nflexcon++;
861             }
862             if (bDynamics || !bFlexCon)
863             {
864                 for (i = 1; i < 3; i++)
865                 {
866                     a = ia[i] - start;
867                     count[a]++;
868                 }
869             }
870             ia += 3;
871         }
872     }
873     *nflexiblecons = nflexcon;
874
875     at2con.nr           = natoms;
876     at2con.nalloc_index = at2con.nr+1;
877     snew(at2con.index, at2con.nalloc_index);
878     at2con.index[0] = 0;
879     for (a = 0; a < natoms; a++)
880     {
881         at2con.index[a+1] = at2con.index[a] + count[a];
882         count[a]          = 0;
883     }
884     at2con.nra      = at2con.index[natoms];
885     at2con.nalloc_a = at2con.nra;
886     snew(at2con.a, at2con.nalloc_a);
887
888     /* The F_CONSTRNC constraints have constraint numbers
889      * that continue after the last F_CONSTR constraint.
890      */
891     con_tot = 0;
892     for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
893     {
894         ncon = ilist[ftype].nr/3;
895         ia   = ilist[ftype].iatoms;
896         for (con = 0; con < ncon; con++)
897         {
898             bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
899                         iparams[ia[0]].constr.dB == 0);
900             if (bDynamics || !bFlexCon)
901             {
902                 for (i = 1; i < 3; i++)
903                 {
904                     a = ia[i] - start;
905                     at2con.a[at2con.index[a]+count[a]++] = con_tot;
906                 }
907             }
908             con_tot++;
909             ia += 3;
910         }
911     }
912
913     sfree(count);
914
915     return at2con;
916 }
917
918 static int *make_at2settle(int natoms, const t_ilist *ilist)
919 {
920     int *at2s;
921     int  a, stride, s;
922
923     snew(at2s, natoms);
924     /* Set all to no settle */
925     for (a = 0; a < natoms; a++)
926     {
927         at2s[a] = -1;
928     }
929
930     stride = 1 + NRAL(F_SETTLE);
931
932     for (s = 0; s < ilist->nr; s += stride)
933     {
934         at2s[ilist->iatoms[s+1]] = s/stride;
935         at2s[ilist->iatoms[s+2]] = s/stride;
936         at2s[ilist->iatoms[s+3]] = s/stride;
937     }
938
939     return at2s;
940 }
941
942 void set_constraints(struct gmx_constr *constr,
943                      gmx_localtop_t *top, t_inputrec *ir,
944                      t_mdatoms *md, t_commrec *cr)
945 {
946     t_idef  *idef;
947     int      ncons;
948     t_ilist *settle;
949     int      iO, iH;
950
951     idef = &top->idef;
952
953     if (constr->ncon_tot > 0)
954     {
955         /* We are using the local topology,
956          * so there are only F_CONSTR constraints.
957          */
958         ncons = idef->il[F_CONSTR].nr/3;
959
960         /* With DD we might also need to call LINCS with ncons=0 for
961          * communicating coordinates to other nodes that do have constraints.
962          */
963         if (ir->eConstrAlg == econtLINCS)
964         {
965             set_lincs(idef, md, EI_DYNAMICS(ir->eI), cr, constr->lincsd);
966         }
967         if (ir->eConstrAlg == econtSHAKE)
968         {
969             if (cr->dd)
970             {
971                 make_shake_sblock_dd(constr, &idef->il[F_CONSTR], &top->cgs, cr->dd);
972             }
973             else
974             {
975                 make_shake_sblock_serial(constr, idef, md);
976             }
977             if (ncons > constr->lagr_nalloc)
978             {
979                 constr->lagr_nalloc = over_alloc_dd(ncons);
980                 srenew(constr->lagr, constr->lagr_nalloc);
981             }
982         }
983     }
984
985     if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
986     {
987         settle          = &idef->il[F_SETTLE];
988         iO              = settle->iatoms[1];
989         iH              = settle->iatoms[2];
990         constr->settled =
991             settle_init(md->massT[iO], md->massT[iH],
992                         md->invmass[iO], md->invmass[iH],
993                         idef->iparams[settle->iatoms[0]].settle.doh,
994                         idef->iparams[settle->iatoms[0]].settle.dhh);
995     }
996
997     /* Make a selection of the local atoms for essential dynamics */
998     if (constr->ed && cr->dd)
999     {
1000         dd_make_local_ed_indices(cr->dd, constr->ed);
1001     }
1002 }
1003
1004 static void constr_recur(t_blocka *at2con,
1005                          t_ilist *ilist, t_iparams *iparams, gmx_bool bTopB,
1006                          int at, int depth, int nc, int *path,
1007                          real r0, real r1, real *r2max,
1008                          int *count)
1009 {
1010     int      ncon1;
1011     t_iatom *ia1, *ia2;
1012     int      c, con, a1;
1013     gmx_bool bUse;
1014     t_iatom *ia;
1015     real     len, rn0, rn1;
1016
1017     (*count)++;
1018
1019     ncon1 = ilist[F_CONSTR].nr/3;
1020     ia1   = ilist[F_CONSTR].iatoms;
1021     ia2   = ilist[F_CONSTRNC].iatoms;
1022
1023     /* Loop over all constraints connected to this atom */
1024     for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
1025     {
1026         con = at2con->a[c];
1027         /* Do not walk over already used constraints */
1028         bUse = TRUE;
1029         for (a1 = 0; a1 < depth; a1++)
1030         {
1031             if (con == path[a1])
1032             {
1033                 bUse = FALSE;
1034             }
1035         }
1036         if (bUse)
1037         {
1038             ia = constr_iatomptr(ncon1, ia1, ia2, con);
1039             /* Flexible constraints currently have length 0, which is incorrect */
1040             if (!bTopB)
1041             {
1042                 len = iparams[ia[0]].constr.dA;
1043             }
1044             else
1045             {
1046                 len = iparams[ia[0]].constr.dB;
1047             }
1048             /* In the worst case the bond directions alternate */
1049             if (nc % 2 == 0)
1050             {
1051                 rn0 = r0 + len;
1052                 rn1 = r1;
1053             }
1054             else
1055             {
1056                 rn0 = r0;
1057                 rn1 = r1 + len;
1058             }
1059             /* Assume angles of 120 degrees between all bonds */
1060             if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
1061             {
1062                 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
1063                 if (debug)
1064                 {
1065                     fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
1066                     for (a1 = 0; a1 < depth; a1++)
1067                     {
1068                         fprintf(debug, " %d %5.3f",
1069                                 path[a1],
1070                                 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
1071                     }
1072                     fprintf(debug, " %d %5.3f\n", con, len);
1073                 }
1074             }
1075             /* Limit the number of recursions to 1000*nc,
1076              * so a call does not take more than a second,
1077              * even for highly connected systems.
1078              */
1079             if (depth + 1 < nc && *count < 1000*nc)
1080             {
1081                 if (ia[1] == at)
1082                 {
1083                     a1 = ia[2];
1084                 }
1085                 else
1086                 {
1087                     a1 = ia[1];
1088                 }
1089                 /* Recursion */
1090                 path[depth] = con;
1091                 constr_recur(at2con, ilist, iparams,
1092                              bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
1093                 path[depth] = -1;
1094             }
1095         }
1096     }
1097 }
1098
1099 static real constr_r_max_moltype(gmx_moltype_t *molt, t_iparams *iparams,
1100                                  t_inputrec *ir)
1101 {
1102     int      natoms, nflexcon, *path, at, count;
1103
1104     t_blocka at2con;
1105     real     r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
1106
1107     if (molt->ilist[F_CONSTR].nr   == 0 &&
1108         molt->ilist[F_CONSTRNC].nr == 0)
1109     {
1110         return 0;
1111     }
1112
1113     natoms = molt->atoms.nr;
1114
1115     at2con = make_at2con(0, natoms, molt->ilist, iparams,
1116                          EI_DYNAMICS(ir->eI), &nflexcon);
1117     snew(path, 1+ir->nProjOrder);
1118     for (at = 0; at < 1+ir->nProjOrder; at++)
1119     {
1120         path[at] = -1;
1121     }
1122
1123     r2maxA = 0;
1124     for (at = 0; at < natoms; at++)
1125     {
1126         r0 = 0;
1127         r1 = 0;
1128
1129         count = 0;
1130         constr_recur(&at2con, molt->ilist, iparams,
1131                      FALSE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxA, &count);
1132     }
1133     if (ir->efep == efepNO)
1134     {
1135         rmax = sqrt(r2maxA);
1136     }
1137     else
1138     {
1139         r2maxB = 0;
1140         for (at = 0; at < natoms; at++)
1141         {
1142             r0    = 0;
1143             r1    = 0;
1144             count = 0;
1145             constr_recur(&at2con, molt->ilist, iparams,
1146                          TRUE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxB, &count);
1147         }
1148         lam0 = ir->fepvals->init_lambda;
1149         if (EI_DYNAMICS(ir->eI))
1150         {
1151             lam0 += ir->init_step*ir->fepvals->delta_lambda;
1152         }
1153         rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
1154         if (EI_DYNAMICS(ir->eI))
1155         {
1156             lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
1157             rmax = max(rmax, (1 - lam1)*sqrt(r2maxA) + lam1*sqrt(r2maxB));
1158         }
1159     }
1160
1161     done_blocka(&at2con);
1162     sfree(path);
1163
1164     return rmax;
1165 }
1166
1167 real constr_r_max(FILE *fplog, gmx_mtop_t *mtop, t_inputrec *ir)
1168 {
1169     int  mt;
1170     real rmax;
1171
1172     rmax = 0;
1173     for (mt = 0; mt < mtop->nmoltype; mt++)
1174     {
1175         rmax = max(rmax,
1176                    constr_r_max_moltype(&mtop->moltype[mt],
1177                                         mtop->ffparams.iparams, ir));
1178     }
1179
1180     if (fplog)
1181     {
1182         fprintf(fplog, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir->nProjOrder, rmax);
1183     }
1184
1185     return rmax;
1186 }
1187
1188 gmx_constr_t init_constraints(FILE *fplog,
1189                               gmx_mtop_t *mtop, t_inputrec *ir,
1190                               gmx_edsam_t ed, t_state *state,
1191                               t_commrec *cr)
1192 {
1193     int                  ncon, nset, nmol, settle_type, i, natoms, mt, nflexcon;
1194     struct gmx_constr   *constr;
1195     char                *env;
1196     t_ilist             *ilist;
1197     gmx_mtop_ilistloop_t iloop;
1198
1199     ncon =
1200         gmx_mtop_ftype_count(mtop, F_CONSTR) +
1201         gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1202     nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
1203
1204     if (ncon+nset == 0 && ir->ePull != epullCONSTRAINT && ed == NULL)
1205     {
1206         return NULL;
1207     }
1208
1209     snew(constr, 1);
1210
1211     constr->ncon_tot = ncon;
1212     constr->nflexcon = 0;
1213     if (ncon > 0)
1214     {
1215         constr->n_at2con_mt = mtop->nmoltype;
1216         snew(constr->at2con_mt, constr->n_at2con_mt);
1217         for (mt = 0; mt < mtop->nmoltype; mt++)
1218         {
1219             constr->at2con_mt[mt] = make_at2con(0, mtop->moltype[mt].atoms.nr,
1220                                                 mtop->moltype[mt].ilist,
1221                                                 mtop->ffparams.iparams,
1222                                                 EI_DYNAMICS(ir->eI), &nflexcon);
1223             for (i = 0; i < mtop->nmolblock; i++)
1224             {
1225                 if (mtop->molblock[i].type == mt)
1226                 {
1227                     constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
1228                 }
1229             }
1230         }
1231
1232         if (constr->nflexcon > 0)
1233         {
1234             if (fplog)
1235             {
1236                 fprintf(fplog, "There are %d flexible constraints\n",
1237                         constr->nflexcon);
1238                 if (ir->fc_stepsize == 0)
1239                 {
1240                     fprintf(fplog, "\n"
1241                             "WARNING: step size for flexible constraining = 0\n"
1242                             "         All flexible constraints will be rigid.\n"
1243                             "         Will try to keep all flexible constraints at their original length,\n"
1244                             "         but the lengths may exhibit some drift.\n\n");
1245                     constr->nflexcon = 0;
1246                 }
1247             }
1248             if (constr->nflexcon > 0)
1249             {
1250                 please_cite(fplog, "Hess2002");
1251             }
1252         }
1253
1254         if (ir->eConstrAlg == econtLINCS)
1255         {
1256             constr->lincsd = init_lincs(fplog, mtop,
1257                                         constr->nflexcon, constr->at2con_mt,
1258                                         DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
1259                                         ir->nLincsIter, ir->nProjOrder);
1260         }
1261
1262         if (ir->eConstrAlg == econtSHAKE)
1263         {
1264             if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1265             {
1266                 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1267             }
1268             if (constr->nflexcon)
1269             {
1270                 gmx_fatal(FARGS, "For this system also velocities and/or forces need to be constrained, this can not be done with SHAKE, you should select LINCS");
1271             }
1272             please_cite(fplog, "Ryckaert77a");
1273             if (ir->bShakeSOR)
1274             {
1275                 please_cite(fplog, "Barth95a");
1276             }
1277
1278             constr->shaked = shake_init();
1279         }
1280     }
1281
1282     if (nset > 0)
1283     {
1284         please_cite(fplog, "Miyamoto92a");
1285
1286         constr->bInterCGsettles = inter_charge_group_settles(mtop);
1287
1288         /* Check that we have only one settle type */
1289         settle_type = -1;
1290         iloop       = gmx_mtop_ilistloop_init(mtop);
1291         while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
1292         {
1293             for (i = 0; i < ilist[F_SETTLE].nr; i += 4)
1294             {
1295                 if (settle_type == -1)
1296                 {
1297                     settle_type = ilist[F_SETTLE].iatoms[i];
1298                 }
1299                 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1300                 {
1301                     gmx_fatal(FARGS,
1302                               "The [molecules] section of your topology specifies more than one block of\n"
1303                               "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1304                               "are trying to partition your solvent into different *groups* (e.g. for\n"
1305                               "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1306                               "files specify groups. Otherwise, you may wish to change the least-used\n"
1307                               "block of molecules with SETTLE constraints into 3 normal constraints.");
1308                 }
1309             }
1310         }
1311
1312         constr->n_at2settle_mt = mtop->nmoltype;
1313         snew(constr->at2settle_mt, constr->n_at2settle_mt);
1314         for (mt = 0; mt < mtop->nmoltype; mt++)
1315         {
1316             constr->at2settle_mt[mt] =
1317                 make_at2settle(mtop->moltype[mt].atoms.nr,
1318                                &mtop->moltype[mt].ilist[F_SETTLE]);
1319         }
1320     }
1321
1322     constr->maxwarn = 999;
1323     env             = getenv("GMX_MAXCONSTRWARN");
1324     if (env)
1325     {
1326         constr->maxwarn = 0;
1327         sscanf(env, "%d", &constr->maxwarn);
1328         if (fplog)
1329         {
1330             fprintf(fplog,
1331                     "Setting the maximum number of constraint warnings to %d\n",
1332                     constr->maxwarn);
1333         }
1334         if (MASTER(cr))
1335         {
1336             fprintf(stderr,
1337                     "Setting the maximum number of constraint warnings to %d\n",
1338                     constr->maxwarn);
1339         }
1340     }
1341     if (constr->maxwarn < 0 && fplog)
1342     {
1343         fprintf(fplog, "maxwarn < 0, will not stop on constraint errors\n");
1344     }
1345     constr->warncount_lincs  = 0;
1346     constr->warncount_settle = 0;
1347
1348     /* Initialize the essential dynamics sampling.
1349      * Put the pointer to the ED struct in constr */
1350     constr->ed = ed;
1351     if (ed != NULL || state->edsamstate.nED > 0)
1352     {
1353         init_edsam(mtop, ir, cr, ed, state->x, state->box, &state->edsamstate);
1354     }
1355
1356     constr->warn_mtop = mtop;
1357
1358     return constr;
1359 }
1360
1361 const t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1362 {
1363     return constr->at2con_mt;
1364 }
1365
1366 const int **atom2settle_moltype(gmx_constr_t constr)
1367 {
1368     return (const int **)constr->at2settle_mt;
1369 }
1370
1371
1372 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop)
1373 {
1374     const gmx_moltype_t *molt;
1375     const t_block       *cgs;
1376     const t_ilist       *il;
1377     int                  mb;
1378     int                  nat, *at2cg, cg, a, ftype, i;
1379     gmx_bool             bInterCG;
1380
1381     bInterCG = FALSE;
1382     for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1383     {
1384         molt = &mtop->moltype[mtop->molblock[mb].type];
1385
1386         if (molt->ilist[F_CONSTR].nr   > 0 ||
1387             molt->ilist[F_CONSTRNC].nr > 0 ||
1388             molt->ilist[F_SETTLE].nr > 0)
1389         {
1390             cgs  = &molt->cgs;
1391             snew(at2cg, molt->atoms.nr);
1392             for (cg = 0; cg < cgs->nr; cg++)
1393             {
1394                 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1395                 {
1396                     at2cg[a] = cg;
1397                 }
1398             }
1399
1400             for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1401             {
1402                 il = &molt->ilist[ftype];
1403                 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(ftype))
1404                 {
1405                     if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1406                     {
1407                         bInterCG = TRUE;
1408                     }
1409                 }
1410             }
1411
1412             sfree(at2cg);
1413         }
1414     }
1415
1416     return bInterCG;
1417 }
1418
1419 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop)
1420 {
1421     const gmx_moltype_t *molt;
1422     const t_block       *cgs;
1423     const t_ilist       *il;
1424     int                  mb;
1425     int                  nat, *at2cg, cg, a, ftype, i;
1426     gmx_bool             bInterCG;
1427
1428     bInterCG = FALSE;
1429     for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1430     {
1431         molt = &mtop->moltype[mtop->molblock[mb].type];
1432
1433         if (molt->ilist[F_SETTLE].nr > 0)
1434         {
1435             cgs  = &molt->cgs;
1436             snew(at2cg, molt->atoms.nr);
1437             for (cg = 0; cg < cgs->nr; cg++)
1438             {
1439                 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1440                 {
1441                     at2cg[a] = cg;
1442                 }
1443             }
1444
1445             for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1446             {
1447                 il = &molt->ilist[ftype];
1448                 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(F_SETTLE))
1449                 {
1450                     if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]] ||
1451                         at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+3]])
1452                     {
1453                         bInterCG = TRUE;
1454                     }
1455                 }
1456             }
1457
1458             sfree(at2cg);
1459         }
1460     }
1461
1462     return bInterCG;
1463 }
1464
1465 /* helper functions for andersen temperature control, because the
1466  * gmx_constr construct is only defined in constr.c. Return the list
1467  * of blocks (get_sblock) and the number of blocks (get_nblocks).  */
1468
1469 extern int *get_sblock(struct gmx_constr *constr)
1470 {
1471     return constr->sblock;
1472 }
1473
1474 extern int get_nblocks(struct gmx_constr *constr)
1475 {
1476     return constr->nblocks;
1477 }