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