0e1a6ee524089e1867cc68db7995b81821178c11
[alexxy/gromacs.git] / src / 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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include "confio.h"
43 #include "constr.h"
44 #include "copyrite.h"
45 #include "invblock.h"
46 #include "main.h"
47 #include "mdrun.h"
48 #include "nrnb.h"
49 #include "smalloc.h"
50 #include "vec.h"
51 #include "physics.h"
52 #include "names.h"
53 #include "txtdump.h"
54 #include "domdec.h"
55 #include "pdbio.h"
56 #include "partdec.h"
57 #include "splitter.h"
58 #include "mtop_util.h"
59 #include "gmxfio.h"
60 #include "gmx_omp_nthreads.h"
61
62 typedef struct gmx_constr {
63     int                ncon_tot;       /* The total number of constraints    */
64     int                nflexcon;       /* The number of flexible constraints */
65     int                n_at2con_mt;    /* The size of at2con = #moltypes     */
66     t_blocka          *at2con_mt;      /* A list of atoms to constraints     */
67     int                n_at2settle_mt; /* The size of at2settle = #moltypes  */
68     int              **at2settle_mt;   /* A list of atoms to settles         */
69     gmx_bool           bInterCGsettles;
70     gmx_lincsdata_t    lincsd;         /* LINCS data                         */
71     gmx_shakedata_t    shaked;         /* SHAKE data                         */
72     gmx_settledata_t   settled;        /* SETTLE data                        */
73     int                nblocks;        /* The number of SHAKE blocks         */
74     int               *sblock;         /* The SHAKE blocks                   */
75     int                sblock_nalloc;  /* The allocation size of sblock      */
76     real              *lagr;           /* Lagrange multipliers for SHAKE     */
77     int                lagr_nalloc;    /* The allocation size of lagr        */
78     int                maxwarn;        /* The maximum number of warnings     */
79     int                warncount_lincs;
80     int                warncount_settle;
81     gmx_edsam_t        ed;            /* The essential dynamics data        */
82
83     tensor            *vir_r_m_dr_th; /* Thread local working data          */
84     int               *settle_error;  /* Thread local working data          */
85
86     gmx_mtop_t        *warn_mtop;     /* Only used for printing warnings    */
87 } t_gmx_constr;
88
89 typedef struct {
90     atom_id iatom[3];
91     atom_id blocknr;
92 } t_sortblock;
93
94 static void *init_vetavars(t_vetavars *vars,
95                            gmx_bool constr_deriv,
96                            real veta, real vetanew, t_inputrec *ir, gmx_ekindata_t *ekind, gmx_bool bPscal)
97 {
98     double g;
99     int    i;
100
101     /* first, set the alpha integrator variable */
102     if ((ir->opts.nrdf[0] > 0) && bPscal)
103     {
104         vars->alpha = 1.0 + DIM/((double)ir->opts.nrdf[0]);
105     }
106     else
107     {
108         vars->alpha = 1.0;
109     }
110     g             = 0.5*veta*ir->delta_t;
111     vars->rscale  = exp(g)*series_sinhx(g);
112     g             = -0.25*vars->alpha*veta*ir->delta_t;
113     vars->vscale  = exp(g)*series_sinhx(g);
114     vars->rvscale = vars->vscale*vars->rscale;
115     vars->veta    = vetanew;
116
117     if (constr_deriv)
118     {
119         snew(vars->vscale_nhc, ir->opts.ngtc);
120         if ((ekind == NULL) || (!bPscal))
121         {
122             for (i = 0; i < ir->opts.ngtc; i++)
123             {
124                 vars->vscale_nhc[i] = 1;
125             }
126         }
127         else
128         {
129             for (i = 0; i < ir->opts.ngtc; i++)
130             {
131                 vars->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
132             }
133         }
134     }
135     else
136     {
137         vars->vscale_nhc = NULL;
138     }
139
140     return vars;
141 }
142
143 static void free_vetavars(t_vetavars *vars)
144 {
145     if (vars->vscale_nhc != NULL)
146     {
147         sfree(vars->vscale_nhc);
148     }
149 }
150
151 static int pcomp(const void *p1, const void *p2)
152 {
153     int          db;
154     atom_id      min1, min2, max1, max2;
155     t_sortblock *a1 = (t_sortblock *)p1;
156     t_sortblock *a2 = (t_sortblock *)p2;
157
158     db = a1->blocknr-a2->blocknr;
159
160     if (db != 0)
161     {
162         return db;
163     }
164
165     min1 = min(a1->iatom[1], a1->iatom[2]);
166     max1 = max(a1->iatom[1], a1->iatom[2]);
167     min2 = min(a2->iatom[1], a2->iatom[2]);
168     max2 = max(a2->iatom[1], a2->iatom[2]);
169
170     if (min1 == min2)
171     {
172         return max1-max2;
173     }
174     else
175     {
176         return min1-min2;
177     }
178 }
179
180 int n_flexible_constraints(struct gmx_constr *constr)
181 {
182     int nflexcon;
183
184     if (constr)
185     {
186         nflexcon = constr->nflexcon;
187     }
188     else
189     {
190         nflexcon = 0;
191     }
192
193     return nflexcon;
194 }
195
196 void too_many_constraint_warnings(int eConstrAlg, int warncount)
197 {
198     const char *abort = "- aborting to avoid logfile runaway.\n"
199         "This normally happens when your system is not sufficiently equilibrated,"
200         "or if you are changing lambda too fast in free energy simulations.\n";
201
202     gmx_fatal(FARGS,
203               "Too many %s warnings (%d)\n"
204               "If you know what you are doing you can %s"
205               "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
206               "but normally it is better to fix the problem",
207               (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
208               (eConstrAlg == econtLINCS) ?
209               "adjust the lincs warning threshold in your mdp file\nor " : "\n");
210 }
211
212 static void write_constr_pdb(const char *fn, const char *title,
213                              gmx_mtop_t *mtop,
214                              int start, int homenr, t_commrec *cr,
215                              rvec x[], matrix box)
216 {
217     char          fname[STRLEN], format[STRLEN];
218     FILE         *out;
219     int           dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
220     gmx_domdec_t *dd;
221     char         *anm, *resnm;
222
223     dd = NULL;
224     if (PAR(cr))
225     {
226         sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
227         if (DOMAINDECOMP(cr))
228         {
229             dd = cr->dd;
230             dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
231             start  = 0;
232             homenr = dd_ac1;
233         }
234     }
235     else
236     {
237         sprintf(fname, "%s.pdb", fn);
238     }
239     sprintf(format, "%s\n", pdbformat);
240
241     out = gmx_fio_fopen(fname, "w");
242
243     fprintf(out, "TITLE     %s\n", title);
244     gmx_write_pdb_box(out, -1, box);
245     for (i = start; i < start+homenr; i++)
246     {
247         if (dd != NULL)
248         {
249             if (i >= dd->nat_home && i < dd_ac0)
250             {
251                 continue;
252             }
253             ii = dd->gatindex[i];
254         }
255         else
256         {
257             ii = i;
258         }
259         gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm);
260         fprintf(out, format, "ATOM", (ii+1)%100000,
261                 anm, resnm, ' ', resnr%10000, ' ',
262                 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ]);
263     }
264     fprintf(out, "TER\n");
265
266     gmx_fio_fclose(out);
267 }
268
269 static void dump_confs(FILE *fplog, gmx_large_int_t step, gmx_mtop_t *mtop,
270                        int start, int homenr, t_commrec *cr,
271                        rvec x[], rvec xprime[], matrix box)
272 {
273     char  buf[256], buf2[22];
274
275     char *env = getenv("GMX_SUPPRESS_DUMP");
276     if (env)
277     {
278         return;
279     }
280
281     sprintf(buf, "step%sb", gmx_step_str(step, buf2));
282     write_constr_pdb(buf, "initial coordinates",
283                      mtop, start, homenr, cr, x, box);
284     sprintf(buf, "step%sc", gmx_step_str(step, buf2));
285     write_constr_pdb(buf, "coordinates after constraining",
286                      mtop, start, homenr, cr, xprime, box);
287     if (fplog)
288     {
289         fprintf(fplog, "Wrote pdb files with previous and current coordinates\n");
290     }
291     fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
292 }
293
294 static void pr_sortblock(FILE *fp, const char *title, int nsb, t_sortblock sb[])
295 {
296     int i;
297
298     fprintf(fp, "%s\n", title);
299     for (i = 0; (i < nsb); i++)
300     {
301         fprintf(fp, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
302                 i, sb[i].iatom[0], sb[i].iatom[1], sb[i].iatom[2],
303                 sb[i].blocknr);
304     }
305 }
306
307 gmx_bool constrain(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
308                    struct gmx_constr *constr,
309                    t_idef *idef, t_inputrec *ir, gmx_ekindata_t *ekind,
310                    t_commrec *cr,
311                    gmx_large_int_t step, int delta_step,
312                    t_mdatoms *md,
313                    rvec *x, rvec *xprime, rvec *min_proj,
314                    gmx_bool bMolPBC, matrix box,
315                    real lambda, real *dvdlambda,
316                    rvec *v, tensor *vir,
317                    t_nrnb *nrnb, int econq, gmx_bool bPscal,
318                    real veta, real vetanew)
319 {
320     gmx_bool    bOK, bDump;
321     int         start, homenr, nrend;
322     int         i, j, d;
323     int         ncons, settle_error;
324     tensor      vir_r_m_dr;
325     rvec       *vstor;
326     real        invdt, vir_fac, t;
327     t_ilist    *settle;
328     int         nsettle;
329     t_pbc       pbc, *pbc_null;
330     char        buf[22];
331     t_vetavars  vetavar;
332     int         nth, th;
333
334     if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
335     {
336         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");
337     }
338
339     bOK   = TRUE;
340     bDump = FALSE;
341
342     start  = md->start;
343     homenr = md->homenr;
344     nrend  = start+homenr;
345
346     /* set constants for pressure control integration */
347     init_vetavars(&vetavar, econq != econqCoord,
348                   veta, vetanew, ir, ekind, bPscal);
349
350     if (ir->delta_t == 0)
351     {
352         invdt = 0;
353     }
354     else
355     {
356         invdt  = 1/ir->delta_t;
357     }
358
359     if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
360     {
361         /* Set the constraint lengths for the step at which this configuration
362          * is meant to be. The invmasses should not be changed.
363          */
364         lambda += delta_step*ir->fepvals->delta_lambda;
365     }
366
367     if (vir != NULL)
368     {
369         clear_mat(vir_r_m_dr);
370     }
371
372     where();
373
374     settle  = &idef->il[F_SETTLE];
375     nsettle = settle->nr/(1+NRAL(F_SETTLE));
376
377     if (nsettle > 0)
378     {
379         nth = gmx_omp_nthreads_get(emntSETTLE);
380     }
381     else
382     {
383         nth = 1;
384     }
385
386     if (nth > 1 && constr->vir_r_m_dr_th == NULL)
387     {
388         snew(constr->vir_r_m_dr_th, nth);
389         snew(constr->settle_error, nth);
390     }
391
392     settle_error = -1;
393
394     /* We do not need full pbc when constraints do not cross charge groups,
395      * i.e. when dd->constraint_comm==NULL.
396      * Note that PBC for constraints is different from PBC for bondeds.
397      * For constraints there is both forward and backward communication.
398      */
399     if (ir->ePBC != epbcNONE &&
400         (cr->dd || bMolPBC) && !(cr->dd && cr->dd->constraint_comm == NULL))
401     {
402         /* With pbc=screw the screw has been changed to a shift
403          * by the constraint coordinate communication routine,
404          * so that here we can use normal pbc.
405          */
406         pbc_null = set_pbc_dd(&pbc, ir->ePBC, cr->dd, FALSE, box);
407     }
408     else
409     {
410         pbc_null = NULL;
411     }
412
413     /* Communicate the coordinates required for the non-local constraints
414      * for LINCS and/or SETTLE.
415      */
416     if (cr->dd)
417     {
418         dd_move_x_constraints(cr->dd, box, x, xprime);
419     }
420     else if (PARTDECOMP(cr))
421     {
422         pd_move_x_constraints(cr, x, xprime);
423     }
424
425     if (constr->lincsd != NULL)
426     {
427         bOK = constrain_lincs(fplog, bLog, bEner, ir, step, constr->lincsd, md, cr,
428                               x, xprime, min_proj,
429                               box, pbc_null, lambda, dvdlambda,
430                               invdt, v, vir != NULL, vir_r_m_dr,
431                               econq, nrnb,
432                               constr->maxwarn, &constr->warncount_lincs);
433         if (!bOK && constr->maxwarn >= 0)
434         {
435             if (fplog != NULL)
436             {
437                 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
438                         econstr_names[econtLINCS], gmx_step_str(step, buf));
439             }
440             bDump = TRUE;
441         }
442     }
443
444     if (constr->nblocks > 0)
445     {
446         switch (econq)
447         {
448             case (econqCoord):
449                 bOK = bshakef(fplog, constr->shaked,
450                               homenr, md->invmass, constr->nblocks, constr->sblock,
451                               idef, ir, x, xprime, nrnb,
452                               constr->lagr, lambda, dvdlambda,
453                               invdt, v, vir != NULL, vir_r_m_dr,
454                               constr->maxwarn >= 0, econq, &vetavar);
455                 break;
456             case (econqVeloc):
457                 bOK = bshakef(fplog, constr->shaked,
458                               homenr, md->invmass, constr->nblocks, constr->sblock,
459                               idef, ir, x, min_proj, nrnb,
460                               constr->lagr, lambda, dvdlambda,
461                               invdt, NULL, vir != NULL, vir_r_m_dr,
462                               constr->maxwarn >= 0, econq, &vetavar);
463                 break;
464             default:
465                 gmx_fatal(FARGS, "Internal error, SHAKE called for constraining something else than coordinates");
466                 break;
467         }
468
469         if (!bOK && constr->maxwarn >= 0)
470         {
471             if (fplog != NULL)
472             {
473                 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
474                         econstr_names[econtSHAKE], gmx_step_str(step, buf));
475             }
476             bDump = TRUE;
477         }
478     }
479
480     if (nsettle > 0)
481     {
482         int calcvir_atom_end;
483
484         if (vir == NULL)
485         {
486             calcvir_atom_end = 0;
487         }
488         else
489         {
490             calcvir_atom_end = md->start + md->homenr;
491         }
492
493         switch (econq)
494         {
495             case econqCoord:
496 #pragma omp parallel for num_threads(nth) schedule(static)
497                 for (th = 0; th < nth; th++)
498                 {
499                     int start_th, end_th;
500
501                     if (th > 0)
502                     {
503                         clear_mat(constr->vir_r_m_dr_th[th]);
504                     }
505
506                     start_th = (nsettle* th   )/nth;
507                     end_th   = (nsettle*(th+1))/nth;
508                     if (start_th >= 0 && end_th - start_th > 0)
509                     {
510                         csettle(constr->settled,
511                                 end_th-start_th,
512                                 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
513                                 pbc_null,
514                                 x[0], xprime[0],
515                                 invdt, v ? v[0] : NULL, calcvir_atom_end,
516                                 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
517                                 th == 0 ? &settle_error : &constr->settle_error[th],
518                                 &vetavar);
519                     }
520                 }
521                 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
522                 if (v != NULL)
523                 {
524                     inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
525                 }
526                 if (vir != NULL)
527                 {
528                     inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
529                 }
530                 break;
531             case econqVeloc:
532             case econqDeriv:
533             case econqForce:
534             case econqForceDispl:
535 #pragma omp parallel for num_threads(nth) schedule(static)
536                 for (th = 0; th < nth; th++)
537                 {
538                     int start_th, end_th;
539
540                     if (th > 0)
541                     {
542                         clear_mat(constr->vir_r_m_dr_th[th]);
543                     }
544
545                     start_th = (nsettle* th   )/nth;
546                     end_th   = (nsettle*(th+1))/nth;
547
548                     if (start_th >= 0 && end_th - start_th > 0)
549                     {
550                         settle_proj(fplog, constr->settled, econq,
551                                     end_th-start_th,
552                                     settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
553                                     pbc_null,
554                                     x,
555                                     xprime, min_proj, calcvir_atom_end,
556                                     th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
557                                     &vetavar);
558                     }
559                 }
560                 /* This is an overestimate */
561                 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
562                 break;
563             case econqDeriv_FlexCon:
564                 /* Nothing to do, since the are no flexible constraints in settles */
565                 break;
566             default:
567                 gmx_incons("Unknown constraint quantity for settle");
568         }
569     }
570
571     if (settle->nr > 0)
572     {
573         /* Combine virial and error info of the other threads */
574         for (i = 1; i < nth; i++)
575         {
576             m_add(vir_r_m_dr, constr->vir_r_m_dr_th[i], vir_r_m_dr);
577             settle_error = constr->settle_error[i];
578         }
579
580         if (econq == econqCoord && settle_error >= 0)
581         {
582             bOK = FALSE;
583             if (constr->maxwarn >= 0)
584             {
585                 char buf[256];
586                 sprintf(buf,
587                         "\nstep " gmx_large_int_pfmt ": Water molecule starting at atom %d can not be "
588                         "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
589                         step, ddglatnr(cr->dd, settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
590                 if (fplog)
591                 {
592                     fprintf(fplog, "%s", buf);
593                 }
594                 fprintf(stderr, "%s", buf);
595                 constr->warncount_settle++;
596                 if (constr->warncount_settle > constr->maxwarn)
597                 {
598                     too_many_constraint_warnings(-1, constr->warncount_settle);
599                 }
600                 bDump = TRUE;
601             }
602         }
603     }
604
605     free_vetavars(&vetavar);
606
607     if (vir != NULL)
608     {
609         switch (econq)
610         {
611             case econqCoord:
612                 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
613                 break;
614             case econqVeloc:
615                 vir_fac = 0.5/ir->delta_t;
616                 break;
617             case econqForce:
618             case econqForceDispl:
619                 vir_fac = 0.5;
620                 break;
621             default:
622                 vir_fac = 0;
623                 gmx_incons("Unsupported constraint quantity for virial");
624         }
625
626         if (EI_VV(ir->eI))
627         {
628             vir_fac *= 2;  /* only constraining over half the distance here */
629         }
630         for (i = 0; i < DIM; i++)
631         {
632             for (j = 0; j < DIM; j++)
633             {
634                 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
635             }
636         }
637     }
638
639     if (bDump)
640     {
641         dump_confs(fplog, step, constr->warn_mtop, start, homenr, cr, x, xprime, box);
642     }
643
644     if (econq == econqCoord)
645     {
646         if (ir->ePull == epullCONSTRAINT)
647         {
648             if (EI_DYNAMICS(ir->eI))
649             {
650                 t = ir->init_t + (step + delta_step)*ir->delta_t;
651             }
652             else
653             {
654                 t = ir->init_t;
655             }
656             set_pbc(&pbc, ir->ePBC, box);
657             pull_constraint(ir->pull, md, &pbc, cr, ir->delta_t, t, x, xprime, v, *vir);
658         }
659         if (constr->ed && delta_step > 0)
660         {
661             /* apply the essential dynamcs constraints here */
662             do_edsam(ir, step, cr, xprime, v, box, constr->ed);
663         }
664     }
665
666     return bOK;
667 }
668
669 real *constr_rmsd_data(struct gmx_constr *constr)
670 {
671     if (constr->lincsd)
672     {
673         return lincs_rmsd_data(constr->lincsd);
674     }
675     else
676     {
677         return NULL;
678     }
679 }
680
681 real constr_rmsd(struct gmx_constr *constr, gmx_bool bSD2)
682 {
683     if (constr->lincsd)
684     {
685         return lincs_rmsd(constr->lincsd, bSD2);
686     }
687     else
688     {
689         return 0;
690     }
691 }
692
693 static void make_shake_sblock_pd(struct gmx_constr *constr,
694                                  t_idef *idef, t_mdatoms *md)
695 {
696     int          i, j, m, ncons;
697     int          bstart, bnr;
698     t_blocka     sblocks;
699     t_sortblock *sb;
700     t_iatom     *iatom;
701     atom_id     *inv_sblock;
702
703     /* Since we are processing the local topology,
704      * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
705      */
706     ncons = idef->il[F_CONSTR].nr/3;
707
708     init_blocka(&sblocks);
709     gen_sblocks(NULL, md->start, md->start+md->homenr, idef, &sblocks, FALSE);
710
711     /*
712        bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
713        nblocks=blocks->multinr[idef->nodeid] - bstart;
714      */
715     bstart          = 0;
716     constr->nblocks = sblocks.nr;
717     if (debug)
718     {
719         fprintf(debug, "ncons: %d, bstart: %d, nblocks: %d\n",
720                 ncons, bstart, constr->nblocks);
721     }
722
723     /* Calculate block number for each atom */
724     inv_sblock = make_invblocka(&sblocks, md->nr);
725
726     done_blocka(&sblocks);
727
728     /* Store the block number in temp array and
729      * sort the constraints in order of the sblock number
730      * and the atom numbers, really sorting a segment of the array!
731      */
732 #ifdef DEBUGIDEF
733     pr_idef(fplog, 0, "Before Sort", idef);
734 #endif
735     iatom = idef->il[F_CONSTR].iatoms;
736     snew(sb, ncons);
737     for (i = 0; (i < ncons); i++, iatom += 3)
738     {
739         for (m = 0; (m < 3); m++)
740         {
741             sb[i].iatom[m] = iatom[m];
742         }
743         sb[i].blocknr = inv_sblock[iatom[1]];
744     }
745
746     /* Now sort the blocks */
747     if (debug)
748     {
749         pr_sortblock(debug, "Before sorting", ncons, sb);
750         fprintf(debug, "Going to sort constraints\n");
751     }
752
753     qsort(sb, ncons, (size_t)sizeof(*sb), pcomp);
754
755     if (debug)
756     {
757         pr_sortblock(debug, "After sorting", ncons, sb);
758     }
759
760     iatom = idef->il[F_CONSTR].iatoms;
761     for (i = 0; (i < ncons); i++, iatom += 3)
762     {
763         for (m = 0; (m < 3); m++)
764         {
765             iatom[m] = sb[i].iatom[m];
766         }
767     }
768 #ifdef DEBUGIDEF
769     pr_idef(fplog, 0, "After Sort", idef);
770 #endif
771
772     j = 0;
773     snew(constr->sblock, constr->nblocks+1);
774     bnr = -2;
775     for (i = 0; (i < ncons); i++)
776     {
777         if (sb[i].blocknr != bnr)
778         {
779             bnr                 = sb[i].blocknr;
780             constr->sblock[j++] = 3*i;
781         }
782     }
783     /* Last block... */
784     constr->sblock[j++] = 3*ncons;
785
786     if (j != (constr->nblocks+1))
787     {
788         fprintf(stderr, "bstart: %d\n", bstart);
789         fprintf(stderr, "j: %d, nblocks: %d, ncons: %d\n",
790                 j, constr->nblocks, ncons);
791         for (i = 0; (i < ncons); i++)
792         {
793             fprintf(stderr, "i: %5d  sb[i].blocknr: %5u\n", i, sb[i].blocknr);
794         }
795         for (j = 0; (j <= constr->nblocks); j++)
796         {
797             fprintf(stderr, "sblock[%3d]=%5d\n", j, (int)constr->sblock[j]);
798         }
799         gmx_fatal(FARGS, "DEATH HORROR: "
800                   "sblocks does not match idef->il[F_CONSTR]");
801     }
802     sfree(sb);
803     sfree(inv_sblock);
804 }
805
806 static void make_shake_sblock_dd(struct gmx_constr *constr,
807                                  t_ilist *ilcon, t_block *cgs,
808                                  gmx_domdec_t *dd)
809 {
810     int      ncons, c, cg;
811     t_iatom *iatom;
812
813     if (dd->ncg_home+1 > constr->sblock_nalloc)
814     {
815         constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
816         srenew(constr->sblock, constr->sblock_nalloc);
817     }
818
819     ncons           = ilcon->nr/3;
820     iatom           = ilcon->iatoms;
821     constr->nblocks = 0;
822     cg              = 0;
823     for (c = 0; c < ncons; c++)
824     {
825         if (c == 0 || iatom[1] >= cgs->index[cg+1])
826         {
827             constr->sblock[constr->nblocks++] = 3*c;
828             while (iatom[1] >= cgs->index[cg+1])
829             {
830                 cg++;
831             }
832         }
833         iatom += 3;
834     }
835     constr->sblock[constr->nblocks] = 3*ncons;
836 }
837
838 t_blocka make_at2con(int start, int natoms,
839                      t_ilist *ilist, t_iparams *iparams,
840                      gmx_bool bDynamics, int *nflexiblecons)
841 {
842     int      *count, ncon, con, con_tot, nflexcon, ftype, i, a;
843     t_iatom  *ia;
844     t_blocka  at2con;
845     gmx_bool  bFlexCon;
846
847     snew(count, natoms);
848     nflexcon = 0;
849     for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
850     {
851         ncon = ilist[ftype].nr/3;
852         ia   = ilist[ftype].iatoms;
853         for (con = 0; con < ncon; con++)
854         {
855             bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
856                         iparams[ia[0]].constr.dB == 0);
857             if (bFlexCon)
858             {
859                 nflexcon++;
860             }
861             if (bDynamics || !bFlexCon)
862             {
863                 for (i = 1; i < 3; i++)
864                 {
865                     a = ia[i] - start;
866                     count[a]++;
867                 }
868             }
869             ia += 3;
870         }
871     }
872     *nflexiblecons = nflexcon;
873
874     at2con.nr           = natoms;
875     at2con.nalloc_index = at2con.nr+1;
876     snew(at2con.index, at2con.nalloc_index);
877     at2con.index[0] = 0;
878     for (a = 0; a < natoms; a++)
879     {
880         at2con.index[a+1] = at2con.index[a] + count[a];
881         count[a]          = 0;
882     }
883     at2con.nra      = at2con.index[natoms];
884     at2con.nalloc_a = at2con.nra;
885     snew(at2con.a, at2con.nalloc_a);
886
887     /* The F_CONSTRNC constraints have constraint numbers
888      * that continue after the last F_CONSTR constraint.
889      */
890     con_tot = 0;
891     for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
892     {
893         ncon = ilist[ftype].nr/3;
894         ia   = ilist[ftype].iatoms;
895         for (con = 0; con < ncon; con++)
896         {
897             bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
898                         iparams[ia[0]].constr.dB == 0);
899             if (bDynamics || !bFlexCon)
900             {
901                 for (i = 1; i < 3; i++)
902                 {
903                     a = ia[i] - start;
904                     at2con.a[at2con.index[a]+count[a]++] = con_tot;
905                 }
906             }
907             con_tot++;
908             ia += 3;
909         }
910     }
911
912     sfree(count);
913
914     return at2con;
915 }
916
917 static int *make_at2settle(int natoms, const t_ilist *ilist)
918 {
919     int *at2s;
920     int  a, stride, s;
921
922     snew(at2s, natoms);
923     /* Set all to no settle */
924     for (a = 0; a < natoms; a++)
925     {
926         at2s[a] = -1;
927     }
928
929     stride = 1 + NRAL(F_SETTLE);
930
931     for (s = 0; s < ilist->nr; s += stride)
932     {
933         at2s[ilist->iatoms[s+1]] = s/stride;
934         at2s[ilist->iatoms[s+2]] = s/stride;
935         at2s[ilist->iatoms[s+3]] = s/stride;
936     }
937
938     return at2s;
939 }
940
941 void set_constraints(struct gmx_constr *constr,
942                      gmx_localtop_t *top, t_inputrec *ir,
943                      t_mdatoms *md, t_commrec *cr)
944 {
945     t_idef  *idef;
946     int      ncons;
947     t_ilist *settle;
948     int      iO, iH;
949
950     idef = &top->idef;
951
952     if (constr->ncon_tot > 0)
953     {
954         /* We are using the local topology,
955          * so there are only F_CONSTR constraints.
956          */
957         ncons = idef->il[F_CONSTR].nr/3;
958
959         /* With DD we might also need to call LINCS with ncons=0 for
960          * communicating coordinates to other nodes that do have constraints.
961          */
962         if (ir->eConstrAlg == econtLINCS)
963         {
964             set_lincs(idef, md, EI_DYNAMICS(ir->eI), cr, constr->lincsd);
965         }
966         if (ir->eConstrAlg == econtSHAKE)
967         {
968             if (cr->dd)
969             {
970                 make_shake_sblock_dd(constr, &idef->il[F_CONSTR], &top->cgs, cr->dd);
971             }
972             else
973             {
974                 make_shake_sblock_pd(constr, idef, md);
975             }
976             if (ncons > constr->lagr_nalloc)
977             {
978                 constr->lagr_nalloc = over_alloc_dd(ncons);
979                 srenew(constr->lagr, constr->lagr_nalloc);
980             }
981         }
982     }
983
984     if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
985     {
986         settle          = &idef->il[F_SETTLE];
987         iO              = settle->iatoms[1];
988         iH              = settle->iatoms[2];
989         constr->settled =
990             settle_init(md->massT[iO], md->massT[iH],
991                         md->invmass[iO], md->invmass[iH],
992                         idef->iparams[settle->iatoms[0]].settle.doh,
993                         idef->iparams[settle->iatoms[0]].settle.dhh);
994     }
995
996     /* Make a selection of the local atoms for essential dynamics */
997     if (constr->ed && cr->dd)
998     {
999         dd_make_local_ed_indices(cr->dd, constr->ed);
1000     }
1001 }
1002
1003 static void constr_recur(t_blocka *at2con,
1004                          t_ilist *ilist, t_iparams *iparams, gmx_bool bTopB,
1005                          int at, int depth, int nc, int *path,
1006                          real r0, real r1, real *r2max,
1007                          int *count)
1008 {
1009     int      ncon1;
1010     t_iatom *ia1, *ia2;
1011     int      c, con, a1;
1012     gmx_bool bUse;
1013     t_iatom *ia;
1014     real     len, rn0, rn1;
1015
1016     (*count)++;
1017
1018     ncon1 = ilist[F_CONSTR].nr/3;
1019     ia1   = ilist[F_CONSTR].iatoms;
1020     ia2   = ilist[F_CONSTRNC].iatoms;
1021
1022     /* Loop over all constraints connected to this atom */
1023     for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
1024     {
1025         con = at2con->a[c];
1026         /* Do not walk over already used constraints */
1027         bUse = TRUE;
1028         for (a1 = 0; a1 < depth; a1++)
1029         {
1030             if (con == path[a1])
1031             {
1032                 bUse = FALSE;
1033             }
1034         }
1035         if (bUse)
1036         {
1037             ia = constr_iatomptr(ncon1, ia1, ia2, con);
1038             /* Flexible constraints currently have length 0, which is incorrect */
1039             if (!bTopB)
1040             {
1041                 len = iparams[ia[0]].constr.dA;
1042             }
1043             else
1044             {
1045                 len = iparams[ia[0]].constr.dB;
1046             }
1047             /* In the worst case the bond directions alternate */
1048             if (nc % 2 == 0)
1049             {
1050                 rn0 = r0 + len;
1051                 rn1 = r1;
1052             }
1053             else
1054             {
1055                 rn0 = r0;
1056                 rn1 = r1 + len;
1057             }
1058             /* Assume angles of 120 degrees between all bonds */
1059             if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
1060             {
1061                 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
1062                 if (debug)
1063                 {
1064                     fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
1065                     for (a1 = 0; a1 < depth; a1++)
1066                     {
1067                         fprintf(debug, " %d %5.3f",
1068                                 path[a1],
1069                                 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
1070                     }
1071                     fprintf(debug, " %d %5.3f\n", con, len);
1072                 }
1073             }
1074             /* Limit the number of recursions to 1000*nc,
1075              * so a call does not take more than a second,
1076              * even for highly connected systems.
1077              */
1078             if (depth + 1 < nc && *count < 1000*nc)
1079             {
1080                 if (ia[1] == at)
1081                 {
1082                     a1 = ia[2];
1083                 }
1084                 else
1085                 {
1086                     a1 = ia[1];
1087                 }
1088                 /* Recursion */
1089                 path[depth] = con;
1090                 constr_recur(at2con, ilist, iparams,
1091                              bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
1092                 path[depth] = -1;
1093             }
1094         }
1095     }
1096 }
1097
1098 static real constr_r_max_moltype(FILE *fplog,
1099                                  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(fplog, &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 }