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