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