Merge remote-tracking branch 'origin/release-4-6' into HEAD
[alexxy/gromacs.git] / src / gromacs / mdlib / shellfc.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11  * Copyright (c) 2001-2008, The GROMACS development team,
12  * check out http://www.gromacs.org for more information.
13  
14  * This program is free software; you can redistribute it and/or
15  * modify it under the terms of the GNU General Public License
16  * as published by the Free Software Foundation; either version 2
17  * of the License, or (at your option) any later version.
18  * 
19  * If you want to redistribute modifications, please consider that
20  * scientific software is very special. Version control is crucial -
21  * bugs must be traceable. We will be happy to consider code for
22  * inclusion in the official distribution, but derived work must not
23  * be called official GROMACS. Details are found in the README & COPYING
24  * files - if they are missing, get the official version at www.gromacs.org.
25  * 
26  * To help us fund GROMACS development, we humbly ask that you cite
27  * the papers on the package - you can find them in the top README file.
28  * 
29  * For more info, check our website at http://www.gromacs.org
30  * 
31  * And Hey:
32  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
33  */
34 #ifdef HAVE_CONFIG_H
35 #include <config.h>
36 #endif
37
38 #include <string.h>
39 #include "typedefs.h"
40 #include "smalloc.h"
41 #include "gmx_fatal.h"
42 #include "vec.h"
43 #include "txtdump.h"
44 #include "mdrun.h"
45 #include "partdec.h"
46 #include "mdatoms.h"
47 #include "vsite.h"
48 #include "network.h"
49 #include "names.h"
50 #include "constr.h"
51 #include "domdec.h"
52 #include "partdec.h"
53 #include "physics.h"
54 #include "copyrite.h"
55 #include "shellfc.h"
56 #include "mtop_util.h"
57 #include "chargegroup.h"
58 #include "macros.h"
59
60
61 typedef struct {
62   int     nnucl;
63   atom_id shell;                /* The shell id                         */
64   atom_id nucl1,nucl2,nucl3;    /* The nuclei connected to the shell    */
65   /* gmx_bool    bInterCG; */       /* Coupled to nuclei outside cg?        */
66   real    k;                    /* force constant                       */
67   real    k_1;                  /* 1 over force constant                */
68   rvec    xold;
69   rvec    fold;
70   rvec    step;
71 } t_shell;
72
73 typedef struct gmx_shellfc {
74   int     nshell_gl;       /* The number of shells in the system       */
75   t_shell *shell_gl;       /* All the shells (for DD only)             */
76   int     *shell_index_gl; /* Global shell index (for DD only)         */
77   gmx_bool    bInterCG;        /* Are there inter charge-group shells?     */
78   int     nshell;          /* The number of local shells               */
79   t_shell *shell;          /* The local shells                         */
80   int     shell_nalloc;    /* The allocation size of shell             */
81   gmx_bool    bPredict;        /* Predict shell positions                  */
82   gmx_bool    bForceInit;      /* Force initialization of shell positions  */
83   int     nflexcon;        /* The number of flexible constraints       */
84   rvec    *x[2];           /* Array for iterative minimization         */
85   rvec    *f[2];           /* Array for iterative minimization         */
86   int     x_nalloc;        /* The allocation size of x and f           */
87   rvec    *acc_dir;        /* Acceleration direction for flexcon       */
88   rvec    *x_old;          /* Old coordinates for flexcon              */
89   int     flex_nalloc;     /* The allocation size of acc_dir and x_old */
90   rvec    *adir_xnold;     /* Work space for init_adir                 */
91   rvec    *adir_xnew;      /* Work space for init_adir                 */
92   int     adir_nalloc;     /* Work space for init_adir                 */
93 } t_gmx_shellfc;
94
95         
96 static void pr_shell(FILE *fplog,int ns,t_shell s[])
97 {
98   int i;
99   
100   fprintf(fplog,"SHELL DATA\n");
101   fprintf(fplog,"%5s  %8s  %5s  %5s  %5s\n",
102           "Shell","Force k","Nucl1","Nucl2","Nucl3");
103   for(i=0; (i<ns); i++) {
104     fprintf(fplog,"%5d  %8.3f  %5d",s[i].shell,1.0/s[i].k_1,s[i].nucl1);
105     if (s[i].nnucl == 2)
106       fprintf(fplog,"  %5d\n",s[i].nucl2);
107     else if (s[i].nnucl == 3)
108       fprintf(fplog,"  %5d  %5d\n",s[i].nucl2,s[i].nucl3);
109     else
110       fprintf(fplog,"\n");
111   }
112 }
113
114 static void predict_shells(FILE *fplog,rvec x[],rvec v[],real dt,
115                            int ns,t_shell s[],
116                            real mass[],gmx_mtop_t *mtop,gmx_bool bInit)
117 {
118   int  i,m,s1,n1,n2,n3;
119   real dt_1,dt_2,dt_3,fudge,tm,m1,m2,m3;
120   rvec *ptr;
121   t_atom *atom;
122   
123   /* We introduce a fudge factor for performance reasons: with this choice
124    * the initial force on the shells is about a factor of two lower than 
125    * without
126    */
127   fudge = 1.0;
128     
129   if (bInit) {
130     if (fplog)
131       fprintf(fplog,"RELAX: Using prediction for initial shell placement\n");
132     ptr  = x;
133     dt_1 = 1;
134   }
135   else {
136     ptr  = v;
137     dt_1 = fudge*dt;
138   }
139     
140   for(i=0; (i<ns); i++) {
141     s1 = s[i].shell;
142     if (bInit)
143       clear_rvec(x[s1]);
144     switch (s[i].nnucl) {
145     case 1:
146       n1 = s[i].nucl1;
147       for(m=0; (m<DIM); m++)
148         x[s1][m]+=ptr[n1][m]*dt_1;
149       break;
150     case 2:
151       n1 = s[i].nucl1;
152       n2 = s[i].nucl2;
153       if (mass) {
154         m1 = mass[n1];
155         m2 = mass[n2];
156       } else {
157         /* Not the correct masses with FE, but it is just a prediction... */
158         m1 = atom[n1].m;
159         m2 = atom[n2].m;
160       }
161       tm = dt_1/(m1+m2);
162       for(m=0; (m<DIM); m++)
163         x[s1][m]+=(m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
164       break;
165     case 3:
166       n1 = s[i].nucl1;
167       n2 = s[i].nucl2;
168       n3 = s[i].nucl3;
169       if (mass) {
170         m1 = mass[n1];
171         m2 = mass[n2];
172         m3 = mass[n3];
173       } else {
174         /* Not the correct masses with FE, but it is just a prediction... */
175         gmx_mtop_atomnr_to_atom(mtop,n1,&atom);
176         m1 = atom->m;
177         gmx_mtop_atomnr_to_atom(mtop,n2,&atom);
178         m2 = atom->m;
179         gmx_mtop_atomnr_to_atom(mtop,n3,&atom);
180         m3 = atom->m;
181       }
182       tm = dt_1/(m1+m2+m3);
183       for(m=0; (m<DIM); m++)
184         x[s1][m]+=(m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
185       break;
186     default:
187       gmx_fatal(FARGS,"Shell %d has %d nuclei!",i,s[i].nnucl);
188     }
189   }
190 }
191
192 gmx_shellfc_t init_shell_flexcon(FILE *fplog,
193                                  gmx_mtop_t *mtop,int nflexcon,
194                                  rvec *x)
195 {
196   struct gmx_shellfc *shfc;
197   t_shell     *shell;
198   int         *shell_index=NULL,*at2cg;
199   t_atom      *atom;
200   int         n[eptNR],ns,nshell,nsi;
201   int         i,j,nmol,type,mb,mt,a_offset,cg,mol,ftype,nra;
202   real        qS,alpha;
203   int         aS,aN=0; /* Shell and nucleus */
204   int         bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
205 #define NBT asize(bondtypes)
206   t_iatom     *ia;
207   gmx_mtop_atomloop_block_t aloopb;
208   gmx_mtop_atomloop_all_t aloop;
209   gmx_ffparams_t *ffparams;
210   gmx_molblock_t *molb;
211   gmx_moltype_t *molt;
212   t_block     *cgs;
213
214   /* Count number of shells, and find their indices */
215   for(i=0; (i<eptNR); i++) {
216     n[i] = 0;
217   }
218
219   aloopb = gmx_mtop_atomloop_block_init(mtop);
220   while (gmx_mtop_atomloop_block_next(aloopb,&atom,&nmol)) {
221     n[atom->ptype] += nmol;
222   }
223
224   if (fplog) {
225     /* Print the number of each particle type */  
226     for(i=0; (i<eptNR); i++) {
227       if (n[i] != 0) {
228         fprintf(fplog,"There are: %d %ss\n",n[i],ptype_str[i]);
229       }
230     }
231   }
232
233   nshell = n[eptShell];
234   
235   if (nshell == 0 && nflexcon == 0) {
236     return NULL;
237   }
238
239   snew(shfc,1);
240   shfc->nflexcon = nflexcon;
241
242   if (nshell == 0) {
243     return shfc;
244   }
245
246   /* We have shells: fill the shell data structure */
247
248   /* Global system sized array, this should be avoided */
249   snew(shell_index,mtop->natoms);
250
251   aloop = gmx_mtop_atomloop_all_init(mtop);
252   nshell = 0;
253   while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
254     if (atom->ptype == eptShell) {
255       shell_index[i] = nshell++;
256     }
257   }
258
259   snew(shell,nshell);
260   
261   /* Initiate the shell structures */    
262   for(i=0; (i<nshell); i++) {
263     shell[i].shell = NO_ATID;
264     shell[i].nnucl = 0;
265     shell[i].nucl1 = NO_ATID;
266     shell[i].nucl2 = NO_ATID;
267     shell[i].nucl3 = NO_ATID;
268     /* shell[i].bInterCG=FALSE; */
269     shell[i].k_1   = 0;
270     shell[i].k     = 0;
271   }
272
273   ffparams = &mtop->ffparams;
274
275   /* Now fill the structures */
276   shfc->bInterCG = FALSE;
277   ns = 0;
278   a_offset = 0;
279   for(mb=0; mb<mtop->nmolblock; mb++) {
280     molb = &mtop->molblock[mb];
281     molt = &mtop->moltype[molb->type];
282
283     cgs = &molt->cgs;
284     snew(at2cg,molt->atoms.nr);
285     for(cg=0; cg<cgs->nr; cg++) {
286       for(i=cgs->index[cg]; i<cgs->index[cg+1]; i++) {
287         at2cg[i] = cg;
288       }
289     }
290
291     atom = molt->atoms.atom;
292     for(mol=0; mol<molb->nmol; mol++) {
293       for(j=0; (j<NBT); j++) {
294         ia = molt->ilist[bondtypes[j]].iatoms;
295         for(i=0; (i<molt->ilist[bondtypes[j]].nr); ) {
296           type  = ia[0];
297           ftype = ffparams->functype[type];
298           nra   = interaction_function[ftype].nratoms;
299           
300           /* Check whether we have a bond with a shell */
301           aS = NO_ATID;
302           
303           switch (bondtypes[j]) {
304           case F_BONDS:
305           case F_HARMONIC:
306           case F_CUBICBONDS:
307           case F_POLARIZATION:
308           case F_ANHARM_POL:
309             if (atom[ia[1]].ptype == eptShell) {
310               aS = ia[1];
311               aN = ia[2];
312             }
313             else if (atom[ia[2]].ptype == eptShell) {
314               aS = ia[2];
315               aN = ia[1];
316             }
317             break;
318           case F_WATER_POL:
319             aN    = ia[4];  /* Dummy */
320             aS    = ia[5];  /* Shell */
321             break;
322           default:
323             gmx_fatal(FARGS,"Death Horror: %s, %d",__FILE__,__LINE__);
324           }
325           
326           if (aS != NO_ATID) {    
327             qS = atom[aS].q;
328             
329             /* Check whether one of the particles is a shell... */
330             nsi = shell_index[a_offset+aS];
331             if ((nsi < 0) || (nsi >= nshell))
332               gmx_fatal(FARGS,"nsi is %d should be within 0 - %d. aS = %d",
333                         nsi,nshell,aS);
334             if (shell[nsi].shell == NO_ATID) {
335               shell[nsi].shell = a_offset + aS;
336               ns ++;
337             }
338             else if (shell[nsi].shell != a_offset+aS)
339               gmx_fatal(FARGS,"Weird stuff in %s, %d",__FILE__,__LINE__);
340             
341             if      (shell[nsi].nucl1 == NO_ATID) {
342               shell[nsi].nucl1 = a_offset + aN;
343             } else if (shell[nsi].nucl2 == NO_ATID) {
344               shell[nsi].nucl2 = a_offset + aN;
345             } else if (shell[nsi].nucl3 == NO_ATID) {
346               shell[nsi].nucl3 = a_offset + aN;
347             } else {
348               if (fplog)
349                 pr_shell(fplog,ns,shell);
350               gmx_fatal(FARGS,"Can not handle more than three bonds per shell\n");
351             }
352             if (at2cg[aS] != at2cg[aN]) {
353               /* shell[nsi].bInterCG = TRUE; */
354               shfc->bInterCG = TRUE;
355             }
356             
357             switch (bondtypes[j]) {
358             case F_BONDS:
359             case F_HARMONIC:
360               shell[nsi].k    += ffparams->iparams[type].harmonic.krA;
361               break;
362             case F_CUBICBONDS:
363               shell[nsi].k    += ffparams->iparams[type].cubic.kb;
364               break;
365             case F_POLARIZATION:
366             case F_ANHARM_POL:
367               if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
368                 gmx_fatal(FARGS,"polarize can not be used with qA(%e) != qB(%e) for atom %d of molecule block %d", qS, atom[aS].qB, aS+1, mb+1);
369               shell[nsi].k    += sqr(qS)*ONE_4PI_EPS0/
370                 ffparams->iparams[type].polarize.alpha;
371               break;
372             case F_WATER_POL:
373               if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
374                 gmx_fatal(FARGS,"water_pol can not be used with qA(%e) != qB(%e) for atom %d of molecule block %d", qS, atom[aS].qB, aS+1, mb+1);
375               alpha          = (ffparams->iparams[type].wpol.al_x+
376                                 ffparams->iparams[type].wpol.al_y+
377                                 ffparams->iparams[type].wpol.al_z)/3.0;
378               shell[nsi].k  += sqr(qS)*ONE_4PI_EPS0/alpha;
379               break;
380             default:
381               gmx_fatal(FARGS,"Death Horror: %s, %d",__FILE__,__LINE__);
382             }
383             shell[nsi].nnucl++;
384           }
385           ia += nra+1;
386           i  += nra+1;
387         }
388       }
389       a_offset += molt->atoms.nr;
390     }
391     /* Done with this molecule type */
392     sfree(at2cg);
393   }
394   
395   /* Verify whether it's all correct */
396   if (ns != nshell)
397     gmx_fatal(FARGS,"Something weird with shells. They may not be bonded to something");
398   
399   for(i=0; (i<ns); i++)
400     shell[i].k_1 = 1.0/shell[i].k;
401   
402   if (debug)
403     pr_shell(debug,ns,shell);
404
405   
406   shfc->nshell_gl      = ns;
407   shfc->shell_gl       = shell;
408   shfc->shell_index_gl = shell_index;
409
410   shfc->bPredict   = (getenv("GMX_NOPREDICT") == NULL);
411   shfc->bForceInit = FALSE;
412   if (!shfc->bPredict) {
413     if (fplog)
414       fprintf(fplog,"\nWill never predict shell positions\n");
415   } else {
416     shfc->bForceInit = (getenv("GMX_FORCEINIT") != NULL);
417     if (shfc->bForceInit && fplog)
418       fprintf(fplog,"\nWill always initiate shell positions\n");
419   }
420
421   if (shfc->bPredict) {
422     if (x) {
423       predict_shells(fplog,x,NULL,0,shfc->nshell_gl,shfc->shell_gl,
424                      NULL,mtop,TRUE);
425     }
426
427     if (shfc->bInterCG) {
428       if (fplog)
429         fprintf(fplog,"\nNOTE: there all shells that are connected to particles outside thier own charge group, will not predict shells positions during the run\n\n");
430       shfc->bPredict = FALSE;
431     }
432   }
433
434   return shfc;
435 }
436
437 void make_local_shells(t_commrec *cr,t_mdatoms *md,
438                        struct gmx_shellfc *shfc)
439 {
440   t_shell *shell;
441   int a0,a1,*ind,nshell,i;
442   gmx_domdec_t *dd=NULL;
443
444   if (PAR(cr)) {
445     if (DOMAINDECOMP(cr)) {
446       dd = cr->dd;
447       a0 = 0;
448       a1 = dd->nat_home;
449     } else {
450       pd_at_range(cr,&a0,&a1);
451     }
452   } else {
453     /* Single node: we need all shells, just copy the pointer */
454     shfc->nshell = shfc->nshell_gl;
455     shfc->shell  = shfc->shell_gl;
456     
457     return;
458   }
459
460   ind = shfc->shell_index_gl;
461
462   nshell = 0;
463   shell  = shfc->shell; 
464   for(i=a0; i<a1; i++) {
465     if (md->ptype[i] == eptShell) {
466       if (nshell+1 > shfc->shell_nalloc) {
467         shfc->shell_nalloc = over_alloc_dd(nshell+1);
468         srenew(shell,shfc->shell_nalloc);
469       }
470       if (dd) {
471         shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
472       } else {
473         shell[nshell] = shfc->shell_gl[ind[i]];
474       }
475       /* With inter-cg shells we can no do shell prediction,
476        * so we do not need the nuclei numbers.
477        */
478       if (!shfc->bInterCG) {
479         shell[nshell].nucl1   = i + shell[nshell].nucl1 - shell[nshell].shell;
480         if (shell[nshell].nnucl > 1)
481           shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
482         if (shell[nshell].nnucl > 2)
483           shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
484       }
485       shell[nshell].shell = i;
486       nshell++;
487     }
488   }
489
490   shfc->nshell = nshell;
491   shfc->shell  = shell;
492 }
493
494 static void do_1pos(rvec xnew,rvec xold,rvec f,real step)
495 {
496   real xo,yo,zo;
497   real dx,dy,dz;
498   
499   xo=xold[XX];
500   yo=xold[YY];
501   zo=xold[ZZ];
502
503   dx=f[XX]*step;
504   dy=f[YY]*step;
505   dz=f[ZZ]*step;
506
507   xnew[XX]=xo+dx;
508   xnew[YY]=yo+dy;
509   xnew[ZZ]=zo+dz;
510 }
511
512 static void do_1pos3(rvec xnew,rvec xold,rvec f,rvec step)
513 {
514   real xo,yo,zo;
515   real dx,dy,dz;
516   
517   xo=xold[XX];
518   yo=xold[YY];
519   zo=xold[ZZ];
520
521   dx=f[XX]*step[XX];
522   dy=f[YY]*step[YY];
523   dz=f[ZZ]*step[ZZ];
524
525   xnew[XX]=xo+dx;
526   xnew[YY]=yo+dy;
527   xnew[ZZ]=zo+dz;
528 }
529
530 static void directional_sd(FILE *log,rvec xold[],rvec xnew[],rvec acc_dir[],
531                            int start,int homenr,real step)
532 {
533   int  i;
534
535   for(i=start; i<homenr; i++)
536     do_1pos(xnew[i],xold[i],acc_dir[i],step);
537 }
538
539 static void shell_pos_sd(FILE *log,rvec xcur[],rvec xnew[],rvec f[],
540                          int ns,t_shell s[],int count)
541 {
542     const real step_scale_min = 0.8,
543         step_scale_increment = 0.2,
544         step_scale_max = 1.2,
545         step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
546   int  i,shell,d;
547   real dx,df,k_est;
548 #ifdef PRINT_STEP  
549   real step_min,step_max;
550
551   step_min = 1e30;
552   step_max = 0;
553 #endif
554   for(i=0; (i<ns); i++) {
555     shell = s[i].shell;
556     if (count == 1) {
557       for(d=0; d<DIM; d++) {
558         s[i].step[d] = s[i].k_1;
559 #ifdef PRINT_STEP
560         step_min = min(step_min,s[i].step[d]);
561         step_max = max(step_max,s[i].step[d]);
562 #endif
563       }
564     } else {
565       for(d=0; d<DIM; d++) {
566         dx = xcur[shell][d] - s[i].xold[d];
567         df =    f[shell][d] - s[i].fold[d];
568     /* -dx/df gets used to generate an interpolated value, but would
569      * cause a NaN if df were binary-equal to zero. Values close to
570      * zero won't cause problems (because of the min() and max()), so
571      * just testing for binary inequality is OK. */
572     if (0.0 != df)
573     {
574         k_est = -dx/df;
575         /* Scale the step size by a factor interpolated from
576          * step_scale_min to step_scale_max, as k_est goes from 0 to
577          * step_scale_multiple * s[i].step[d] */
578         s[i].step[d] =
579             step_scale_min * s[i].step[d] +
580             step_scale_increment * min(step_scale_multiple * s[i].step[d], max(k_est, 0));
581     }
582     else
583     {
584         /* Here 0 == df */
585         if (gmx_numzero(dx)) /* 0 == dx */
586         {
587             /* Likely this will never happen, but if it does just
588              * don't scale the step. */
589         }
590         else /* 0 != dx */
591         {
592             s[i].step[d] *= step_scale_max;
593         }
594     }
595 #ifdef PRINT_STEP
596         step_min = min(step_min,s[i].step[d]);
597         step_max = max(step_max,s[i].step[d]);
598 #endif
599       }
600     }
601     copy_rvec(xcur[shell],s[i].xold);
602     copy_rvec(f[shell],   s[i].fold);
603
604     do_1pos3(xnew[shell],xcur[shell],f[shell],s[i].step);
605
606     if (gmx_debug_at) {
607       fprintf(debug,"shell[%d] = %d\n",i,shell);
608       pr_rvec(debug,0,"fshell",f[shell],DIM,TRUE);
609       pr_rvec(debug,0,"xold",xcur[shell],DIM,TRUE);
610       pr_rvec(debug,0,"step",s[i].step,DIM,TRUE);
611       pr_rvec(debug,0,"xnew",xnew[shell],DIM,TRUE);
612     }
613   }
614 #ifdef PRINT_STEP
615   printf("step %.3e %.3e\n",step_min,step_max);
616 #endif
617 }
618
619 static void decrease_step_size(int nshell,t_shell s[])
620 {
621   int i;
622   
623   for(i=0; i<nshell; i++)
624     svmul(0.8,s[i].step,s[i].step);
625 }
626
627 static void print_epot(FILE *fp,gmx_large_int_t mdstep,int count,real epot,real df,
628                        int ndir,real sf_dir)
629 {
630   char buf[22];
631
632   fprintf(fp,"MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
633           gmx_step_str(mdstep,buf),count,epot,df);
634   if (ndir)
635     fprintf(fp,", dir. rmsF: %6.2e\n",sqrt(sf_dir/ndir));
636   else
637     fprintf(fp,"\n");
638 }
639
640
641 static real rms_force(t_commrec *cr,rvec f[],int ns,t_shell s[],
642                       int ndir,real *sf_dir,real *Epot)
643 {
644   int  i,shell,ntot;
645   double buf[4];
646
647   buf[0] = *sf_dir;
648   for(i=0; i<ns; i++) {
649     shell = s[i].shell;
650     buf[0]  += norm2(f[shell]);
651   }
652   ntot = ns;
653
654   if (PAR(cr)) {
655     buf[1] = ntot;
656     buf[2] = *sf_dir;
657     buf[3] = *Epot;
658     gmx_sumd(4,buf,cr);
659     ntot = (int)(buf[1] + 0.5);
660     *sf_dir = buf[2];
661     *Epot   = buf[3];
662   }
663   ntot += ndir;
664
665   return (ntot ? sqrt(buf[0]/ntot) : 0);
666 }
667
668 static void check_pbc(FILE *fp,rvec x[],int shell)
669 {
670   int m,now;
671   
672   now = shell-4;
673   for(m=0; (m<DIM); m++)
674     if (fabs(x[shell][m]-x[now][m]) > 0.3) {
675       pr_rvecs(fp,0,"SHELL-X",x+now,5);
676       break;
677     }
678 }
679
680 static void dump_shells(FILE *fp,rvec x[],rvec f[],real ftol,int ns,t_shell s[])
681 {
682   int  i,shell;
683   real ft2,ff2;
684   
685   ft2 = sqr(ftol);
686   
687   for(i=0; (i<ns); i++) {
688     shell = s[i].shell;
689     ff2   = iprod(f[shell],f[shell]);
690     if (ff2 > ft2)
691       fprintf(fp,"SHELL %5d, force %10.5f  %10.5f  %10.5f, |f| %10.5f\n",
692               shell,f[shell][XX],f[shell][YY],f[shell][ZZ],sqrt(ff2));
693     check_pbc(fp,x,shell);
694   }
695 }
696
697 static void init_adir(FILE *log,gmx_shellfc_t shfc,
698                       gmx_constr_t constr,t_idef *idef,t_inputrec *ir,
699                       t_commrec *cr,int dd_ac1,
700                       gmx_large_int_t step,t_mdatoms *md,int start,int end,
701                       rvec *x_old,rvec *x_init,rvec *x,
702                       rvec *f,rvec *acc_dir,matrix box,
703                       real *lambda,real *dvdlambda,t_nrnb *nrnb)
704 {
705   rvec   *xnold,*xnew;
706   double w_dt;
707   int    gf,ga,gt;
708   real   dt,scale;
709   int    n,d; 
710   unsigned short *ptype;
711   rvec   p,dx;
712   
713   if (DOMAINDECOMP(cr))
714     n = dd_ac1;
715   else
716     n = end - start;
717   if (n > shfc->adir_nalloc) {
718     shfc->adir_nalloc = over_alloc_dd(n);
719     srenew(shfc->adir_xnold,shfc->adir_nalloc);
720     srenew(shfc->adir_xnew ,shfc->adir_nalloc);
721   }
722   xnold = shfc->adir_xnold;
723   xnew  = shfc->adir_xnew;
724     
725   ptype = md->ptype;
726
727   dt = ir->delta_t;
728
729   /* Does NOT work with freeze or acceleration groups (yet) */
730   for (n=start; n<end; n++) {  
731     w_dt = md->invmass[n]*dt;
732     
733     for (d=0; d<DIM; d++) {
734       if ((ptype[n] != eptVSite) && (ptype[n] != eptShell)) {
735         xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
736         xnew[n-start][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
737       } else {
738         xnold[n-start][d] = x[n][d];
739         xnew[n-start][d] = x[n][d];
740       }
741     }
742   }
743   constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
744             x,xnold-start,NULL,box,
745             lambda[efptBONDED],&(dvdlambda[efptBONDED]),NULL,NULL,nrnb,econqCoord,FALSE,0,0);
746   constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
747             x,xnew-start,NULL,box,
748             lambda[efptBONDED],&(dvdlambda[efptBONDED]),NULL,NULL,nrnb,econqCoord,FALSE,0,0);
749
750   /* Set xnew to minus the acceleration */
751   for (n=start; n<end; n++) {
752     for(d=0; d<DIM; d++)
753       xnew[n-start][d] =
754         -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
755         - f[n][d]*md->invmass[n];
756     clear_rvec(acc_dir[n]);
757   }
758
759   /* Project the acceleration on the old bond directions */
760   constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
761             x_old,xnew-start,acc_dir,box,
762             lambda[efptBONDED],&(dvdlambda[efptBONDED]),NULL,NULL,nrnb,econqDeriv_FlexCon,FALSE,0,0);
763 }
764
765 int relax_shell_flexcon(FILE *fplog,t_commrec *cr,gmx_bool bVerbose,
766                         gmx_large_int_t mdstep,t_inputrec *inputrec,
767                         gmx_bool bDoNS,int force_flags,
768                         gmx_bool bStopCM,
769                         gmx_localtop_t *top,
770                         gmx_mtop_t* mtop,
771                         gmx_constr_t constr,
772                         gmx_enerdata_t *enerd,t_fcdata *fcd,
773                         t_state *state,rvec f[],
774                         tensor force_vir,
775                         t_mdatoms *md,
776                         t_nrnb *nrnb,gmx_wallcycle_t wcycle,
777                         t_graph *graph,
778                         gmx_groups_t *groups,
779                         struct gmx_shellfc *shfc,
780                         t_forcerec *fr,
781                         gmx_bool bBornRadii,
782                         double t,rvec mu_tot,
783                         int natoms,gmx_bool *bConverged,
784                         gmx_vsite_t *vsite,
785                         FILE *fp_field)
786 {
787   int    nshell;
788   t_shell *shell;
789   t_idef *idef;
790   rvec   *pos[2],*force[2],*acc_dir=NULL,*x_old=NULL;
791   real   Epot[2],df[2];
792   rvec   dx;
793   real   sf_dir,invdt;
794   real   ftol,xiH,xiS,dum=0;
795   char   sbuf[22];
796   gmx_bool   bCont,bInit;
797   int    nat,dd_ac0,dd_ac1=0,i;
798   int    start=md->start,homenr=md->homenr,end=start+homenr,cg0,cg1;
799   int    nflexcon,g,number_steps,d,Min=0,count=0;
800 #define  Try (1-Min)             /* At start Try = 1 */
801
802   bCont        = (mdstep == inputrec->init_step) && inputrec->bContinuation;
803   bInit        = (mdstep == inputrec->init_step) || shfc->bForceInit;
804   ftol         = inputrec->em_tol;
805   number_steps = inputrec->niter;
806   nshell       = shfc->nshell;
807   shell        = shfc->shell;
808   nflexcon     = shfc->nflexcon;
809
810   idef = &top->idef;
811
812   if (DOMAINDECOMP(cr)) {
813     nat = dd_natoms_vsite(cr->dd);
814     if (nflexcon > 0) {
815       dd_get_constraint_range(cr->dd,&dd_ac0,&dd_ac1);
816       nat = max(nat,dd_ac1);
817     }
818   } else {
819     nat = state->natoms;
820   }
821
822   if (nat > shfc->x_nalloc) {
823     /* Allocate local arrays */
824     shfc->x_nalloc = over_alloc_dd(nat);
825     for(i=0; (i<2); i++) {
826       srenew(shfc->x[i],shfc->x_nalloc);
827       srenew(shfc->f[i],shfc->x_nalloc);
828     }
829   }
830   for(i=0; (i<2); i++) {
831     pos[i]   = shfc->x[i];
832     force[i] = shfc->f[i];
833   }
834      
835   /* With particle decomposition this code only works
836    * when all particles involved with each shell are in the same cg.
837    */
838
839   if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr)) {
840     /* This is the only time where the coordinates are used
841      * before do_force is called, which normally puts all
842      * charge groups in the box.
843      */
844     if (PARTDECOMP(cr)) {
845       pd_cg_range(cr,&cg0,&cg1);
846     } else {
847       cg0 = 0;
848       cg1 = top->cgs.nr;
849     }
850     put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,state->box,
851                              &(top->cgs),state->x,fr->cg_cm);
852     if (graph)
853       mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
854   }
855
856   /* After this all coordinate arrays will contain whole molecules */
857   if (graph)
858     shift_self(graph,state->box,state->x);
859
860   if (nflexcon) {
861     if (nat > shfc->flex_nalloc) {
862       shfc->flex_nalloc = over_alloc_dd(nat);
863       srenew(shfc->acc_dir,shfc->flex_nalloc);
864       srenew(shfc->x_old,shfc->flex_nalloc);
865     }
866     acc_dir = shfc->acc_dir;
867     x_old   = shfc->x_old;
868     for(i=0; i<homenr; i++) {
869       for(d=0; d<DIM; d++)
870         shfc->x_old[i][d] =
871           state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
872     }
873   }
874
875   /* Do a prediction of the shell positions */
876   if (shfc->bPredict && !bCont) {
877     predict_shells(fplog,state->x,state->v,inputrec->delta_t,nshell,shell,
878                    md->massT,NULL,bInit);
879   }
880
881   /* do_force expected the charge groups to be in the box */
882   if (graph)
883     unshift_self(graph,state->box,state->x);
884
885   /* Calculate the forces first time around */
886   if (gmx_debug_at) {
887     pr_rvecs(debug,0,"x b4 do_force",state->x + start,homenr);
888   }
889   do_force(fplog,cr,inputrec,mdstep,nrnb,wcycle,top,mtop,groups,
890            state->box,state->x,&state->hist,
891            force[Min],force_vir,md,enerd,fcd,
892            state->lambda,graph,
893            fr,vsite,mu_tot,t,fp_field,NULL,bBornRadii,
894            (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
895
896   sf_dir = 0;
897   if (nflexcon) {
898     init_adir(fplog,shfc,
899               constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
900               shfc->x_old-start,state->x,state->x,force[Min],
901               shfc->acc_dir-start,state->box,state->lambda,&dum,nrnb);
902
903     for(i=start; i<end; i++)
904       sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
905   }
906
907   Epot[Min] = enerd->term[F_EPOT];
908
909   df[Min]=rms_force(cr,shfc->f[Min],nshell,shell,nflexcon,&sf_dir,&Epot[Min]);
910   df[Try]=0;
911   if (debug) {
912     fprintf(debug,"df = %g  %g\n",df[Min],df[Try]);
913   }
914
915   if (gmx_debug_at) {
916     pr_rvecs(debug,0,"force0",force[Min],md->nr);
917   }
918
919   if (nshell+nflexcon > 0) {
920     /* Copy x to pos[Min] & pos[Try]: during minimization only the
921      * shell positions are updated, therefore the other particles must
922      * be set here.
923      */
924     memcpy(pos[Min],state->x,nat*sizeof(state->x[0]));
925     memcpy(pos[Try],state->x,nat*sizeof(state->x[0]));
926   }
927   
928   if (bVerbose && MASTER(cr))
929     print_epot(stdout,mdstep,0,Epot[Min],df[Min],nflexcon,sf_dir);
930
931   if (debug) {
932     fprintf(debug,"%17s: %14.10e\n",
933             interaction_function[F_EKIN].longname,enerd->term[F_EKIN]);
934     fprintf(debug,"%17s: %14.10e\n",
935             interaction_function[F_EPOT].longname,enerd->term[F_EPOT]);
936     fprintf(debug,"%17s: %14.10e\n",
937             interaction_function[F_ETOT].longname,enerd->term[F_ETOT]);
938     fprintf(debug,"SHELLSTEP %s\n",gmx_step_str(mdstep,sbuf));
939   }
940   
941   /* First check whether we should do shells, or whether the force is 
942    * low enough even without minimization.
943    */
944   *bConverged = (df[Min] < ftol);
945   
946   for(count=1; (!(*bConverged) && (count < number_steps)); count++) {
947     if (vsite)
948       construct_vsites(fplog,vsite,pos[Min],nrnb,inputrec->delta_t,state->v,
949                        idef->iparams,idef->il,
950                        fr->ePBC,fr->bMolPBC,graph,cr,state->box);
951      
952     if (nflexcon) {
953       init_adir(fplog,shfc,
954                 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
955                 x_old-start,state->x,pos[Min],force[Min],acc_dir-start,
956                 state->box,state->lambda,&dum,nrnb);
957       
958       directional_sd(fplog,pos[Min],pos[Try],acc_dir-start,start,end,
959                      fr->fc_stepsize);
960     }
961     
962     /* New positions, Steepest descent */
963     shell_pos_sd(fplog,pos[Min],pos[Try],force[Min],nshell,shell,count); 
964
965     /* do_force expected the charge groups to be in the box */
966     if (graph)
967       unshift_self(graph,state->box,pos[Try]);
968
969     if (gmx_debug_at) {
970       pr_rvecs(debug,0,"RELAX: pos[Min]  ",pos[Min] + start,homenr);
971       pr_rvecs(debug,0,"RELAX: pos[Try]  ",pos[Try] + start,homenr);
972     }
973     /* Try the new positions */
974     do_force(fplog,cr,inputrec,1,nrnb,wcycle,
975              top,mtop,groups,state->box,pos[Try],&state->hist,
976              force[Try],force_vir,
977              md,enerd,fcd,state->lambda,graph,
978              fr,vsite,mu_tot,t,fp_field,NULL,bBornRadii,
979              force_flags);
980     
981     if (gmx_debug_at) {
982       pr_rvecs(debug,0,"RELAX: force[Min]",force[Min] + start,homenr);
983       pr_rvecs(debug,0,"RELAX: force[Try]",force[Try] + start,homenr);
984     }
985     sf_dir = 0;
986     if (nflexcon) {
987       init_adir(fplog,shfc,
988                 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
989                 x_old-start,state->x,pos[Try],force[Try],acc_dir-start,
990                 state->box,state->lambda,&dum,nrnb);
991
992       for(i=start; i<end; i++)
993         sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
994     }
995
996     Epot[Try] = enerd->term[F_EPOT]; 
997     
998     df[Try]=rms_force(cr,force[Try],nshell,shell,nflexcon,&sf_dir,&Epot[Try]);
999
1000     if (debug)
1001       fprintf(debug,"df = %g  %g\n",df[Min],df[Try]);
1002
1003     if (debug) {
1004       if (gmx_debug_at)
1005         pr_rvecs(debug,0,"F na do_force",force[Try] + start,homenr);
1006       if (gmx_debug_at) {
1007         fprintf(debug,"SHELL ITER %d\n",count);
1008         dump_shells(debug,pos[Try],force[Try],ftol,nshell,shell);
1009       }
1010     }
1011
1012     if (bVerbose && MASTER(cr))
1013       print_epot(stdout,mdstep,count,Epot[Try],df[Try],nflexcon,sf_dir);
1014       
1015     *bConverged = (df[Try] < ftol);
1016     
1017     if ((df[Try] < df[Min])) {
1018       if (debug)
1019         fprintf(debug,"Swapping Min and Try\n");
1020       if (nflexcon) {
1021         /* Correct the velocities for the flexible constraints */
1022         invdt = 1/inputrec->delta_t;
1023         for(i=start; i<end; i++) {
1024           for(d=0; d<DIM; d++)
1025             state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1026         }
1027       }
1028       Min  = Try;
1029     } else {
1030       decrease_step_size(nshell,shell);
1031     }
1032   }
1033   if (MASTER(cr) && !(*bConverged)) {
1034     /* Note that the energies and virial are incorrect when not converged */
1035     if (fplog)
1036       fprintf(fplog,
1037               "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1038               gmx_step_str(mdstep,sbuf),number_steps,df[Min]);
1039     fprintf(stderr,
1040             "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1041             gmx_step_str(mdstep,sbuf),number_steps,df[Min]);
1042   }
1043
1044   /* Copy back the coordinates and the forces */
1045   memcpy(state->x,pos[Min],nat*sizeof(state->x[0]));
1046   memcpy(f,force[Min],nat*sizeof(f[0]));
1047
1048   return count; 
1049 }
1050