Added two new bonded functions.
[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_ANHARM_POL, 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           case F_ANHARM_POL:
308             if (atom[ia[1]].ptype == eptShell) {
309               aS = ia[1];
310               aN = ia[2];
311             }
312             else if (atom[ia[2]].ptype == eptShell) {
313               aS = ia[2];
314               aN = ia[1];
315             }
316             break;
317           case F_WATER_POL:
318             aN    = ia[4];  /* Dummy */
319             aS    = ia[5];  /* Shell */
320             break;
321           default:
322             gmx_fatal(FARGS,"Death Horror: %s, %d",__FILE__,__LINE__);
323           }
324           
325           if (aS != NO_ATID) {    
326             qS = atom[aS].q;
327             
328             /* Check whether one of the particles is a shell... */
329             nsi = shell_index[a_offset+aS];
330             if ((nsi < 0) || (nsi >= nshell))
331               gmx_fatal(FARGS,"nsi is %d should be within 0 - %d. aS = %d",
332                         nsi,nshell,aS);
333             if (shell[nsi].shell == NO_ATID) {
334               shell[nsi].shell = a_offset + aS;
335               ns ++;
336             }
337             else if (shell[nsi].shell != a_offset+aS)
338               gmx_fatal(FARGS,"Weird stuff in %s, %d",__FILE__,__LINE__);
339             
340             if      (shell[nsi].nucl1 == NO_ATID) {
341               shell[nsi].nucl1 = a_offset + aN;
342             } else if (shell[nsi].nucl2 == NO_ATID) {
343               shell[nsi].nucl2 = a_offset + aN;
344             } else if (shell[nsi].nucl3 == NO_ATID) {
345               shell[nsi].nucl3 = a_offset + aN;
346             } else {
347               if (fplog)
348                 pr_shell(fplog,ns,shell);
349               gmx_fatal(FARGS,"Can not handle more than three bonds per shell\n");
350             }
351             if (at2cg[aS] != at2cg[aN]) {
352               /* shell[nsi].bInterCG = TRUE; */
353               shfc->bInterCG = TRUE;
354             }
355             
356             switch (bondtypes[j]) {
357             case F_BONDS:
358             case F_HARMONIC:
359               shell[nsi].k    += ffparams->iparams[type].harmonic.krA;
360               break;
361             case F_CUBICBONDS:
362               shell[nsi].k    += ffparams->iparams[type].cubic.kb;
363               break;
364             case F_POLARIZATION:
365             case F_ANHARM_POL:
366               if (qS != atom[aS].qB)
367                 gmx_fatal(FARGS,"polarize can not be used with qA != qB");
368               shell[nsi].k    += sqr(qS)*ONE_4PI_EPS0/
369                 ffparams->iparams[type].polarize.alpha;
370               break;
371             case F_WATER_POL:
372               if (qS != atom[aS].qB)
373                 gmx_fatal(FARGS,"water_pol can not be used with qA != qB");
374               alpha          = (ffparams->iparams[type].wpol.al_x+
375                                 ffparams->iparams[type].wpol.al_y+
376                                 ffparams->iparams[type].wpol.al_z)/3.0;
377               shell[nsi].k  += sqr(qS)*ONE_4PI_EPS0/alpha;
378               break;
379             default:
380               gmx_fatal(FARGS,"Death Horror: %s, %d",__FILE__,__LINE__);
381             }
382             shell[nsi].nnucl++;
383           }
384           ia += nra+1;
385           i  += nra+1;
386         }
387       }
388       a_offset += molt->atoms.nr;
389     }
390     /* Done with this molecule type */
391     sfree(at2cg);
392   }
393   
394   /* Verify whether it's all correct */
395   if (ns != nshell)
396     gmx_fatal(FARGS,"Something weird with shells. They may not be bonded to something");
397   
398   for(i=0; (i<ns); i++)
399     shell[i].k_1 = 1.0/shell[i].k;
400   
401   if (debug)
402     pr_shell(debug,ns,shell);
403
404   
405   shfc->nshell_gl      = ns;
406   shfc->shell_gl       = shell;
407   shfc->shell_index_gl = shell_index;
408
409   shfc->bPredict   = (getenv("GMX_NOPREDICT") == NULL);
410   shfc->bForceInit = FALSE;
411   if (!shfc->bPredict) {
412     if (fplog)
413       fprintf(fplog,"\nWill never predict shell positions\n");
414   } else {
415     shfc->bForceInit = (getenv("GMX_FORCEINIT") != NULL);
416     if (shfc->bForceInit && fplog)
417       fprintf(fplog,"\nWill always initiate shell positions\n");
418   }
419
420   if (shfc->bPredict) {
421     if (x) {
422       predict_shells(fplog,x,NULL,0,shfc->nshell_gl,shfc->shell_gl,
423                      NULL,mtop,TRUE);
424     }
425
426     if (shfc->bInterCG) {
427       if (fplog)
428         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");
429       shfc->bPredict = FALSE;
430     }
431   }
432
433   return shfc;
434 }
435
436 void make_local_shells(t_commrec *cr,t_mdatoms *md,
437                        struct gmx_shellfc *shfc)
438 {
439   t_shell *shell;
440   int a0,a1,*ind,nshell,i;
441   gmx_domdec_t *dd=NULL;
442
443   if (PAR(cr)) {
444     if (DOMAINDECOMP(cr)) {
445       dd = cr->dd;
446       a0 = 0;
447       a1 = dd->nat_home;
448     } else {
449       pd_at_range(cr,&a0,&a1);
450     }
451   } else {
452     /* Single node: we need all shells, just copy the pointer */
453     shfc->nshell = shfc->nshell_gl;
454     shfc->shell  = shfc->shell_gl;
455     
456     return;
457   }
458
459   ind = shfc->shell_index_gl;
460
461   nshell = 0;
462   shell  = shfc->shell; 
463   for(i=a0; i<a1; i++) {
464     if (md->ptype[i] == eptShell) {
465       if (nshell+1 > shfc->shell_nalloc) {
466         shfc->shell_nalloc = over_alloc_dd(nshell+1);
467         srenew(shell,shfc->shell_nalloc);
468       }
469       if (dd) {
470         shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
471       } else {
472         shell[nshell] = shfc->shell_gl[ind[i]];
473       }
474       /* With inter-cg shells we can no do shell prediction,
475        * so we do not need the nuclei numbers.
476        */
477       if (!shfc->bInterCG) {
478         shell[nshell].nucl1   = i + shell[nshell].nucl1 - shell[nshell].shell;
479         if (shell[nshell].nnucl > 1)
480           shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
481         if (shell[nshell].nnucl > 2)
482           shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
483       }
484       shell[nshell].shell = i;
485       nshell++;
486     }
487   }
488
489   shfc->nshell = nshell;
490   shfc->shell  = shell;
491 }
492
493 static void do_1pos(rvec xnew,rvec xold,rvec f,real step)
494 {
495   real xo,yo,zo;
496   real dx,dy,dz;
497   
498   xo=xold[XX];
499   yo=xold[YY];
500   zo=xold[ZZ];
501
502   dx=f[XX]*step;
503   dy=f[YY]*step;
504   dz=f[ZZ]*step;
505
506   xnew[XX]=xo+dx;
507   xnew[YY]=yo+dy;
508   xnew[ZZ]=zo+dz;
509 }
510
511 static void do_1pos3(rvec xnew,rvec xold,rvec f,rvec step)
512 {
513   real xo,yo,zo;
514   real dx,dy,dz;
515   
516   xo=xold[XX];
517   yo=xold[YY];
518   zo=xold[ZZ];
519
520   dx=f[XX]*step[XX];
521   dy=f[YY]*step[YY];
522   dz=f[ZZ]*step[ZZ];
523
524   xnew[XX]=xo+dx;
525   xnew[YY]=yo+dy;
526   xnew[ZZ]=zo+dz;
527 }
528
529 static void directional_sd(FILE *log,rvec xold[],rvec xnew[],rvec acc_dir[],
530                            int start,int homenr,real step)
531 {
532   int  i;
533
534   for(i=start; i<homenr; i++)
535     do_1pos(xnew[i],xold[i],acc_dir[i],step);
536 }
537
538 static void shell_pos_sd(FILE *log,rvec xcur[],rvec xnew[],rvec f[],
539                          int ns,t_shell s[],int count)
540 {
541   int  i,shell,d;
542   real dx,df,k_est;
543 #ifdef PRINT_STEP  
544   real step_min,step_max;
545
546   step_min = 1e30;
547   step_max = 0;
548 #endif
549   for(i=0; (i<ns); i++) {
550     shell = s[i].shell;
551     if (count == 1) {
552       for(d=0; d<DIM; d++) {
553         s[i].step[d] = s[i].k_1;
554 #ifdef PRINT_STEP
555         step_min = min(step_min,s[i].step[d]);
556         step_max = max(step_max,s[i].step[d]);
557 #endif
558       }
559     } else {
560       for(d=0; d<DIM; d++) {
561         dx = xcur[shell][d] - s[i].xold[d];
562         df =    f[shell][d] - s[i].fold[d];
563         if (dx != 0 && df != 0) {
564           k_est = -dx/df;
565           if (k_est >= 2*s[i].step[d]) {
566             s[i].step[d] *= 1.2;
567           } else if (k_est <= 0) {
568             s[i].step[d] *= 0.8;
569           } else {
570             s[i].step[d] = 0.8*s[i].step[d] + 0.2*k_est;
571           }
572         } else if (dx != 0) {
573           s[i].step[d] *= 1.2;
574         }
575 #ifdef PRINT_STEP
576         step_min = min(step_min,s[i].step[d]);
577         step_max = max(step_max,s[i].step[d]);
578 #endif
579       }
580     }
581     copy_rvec(xcur[shell],s[i].xold);
582     copy_rvec(f[shell],   s[i].fold);
583
584     do_1pos3(xnew[shell],xcur[shell],f[shell],s[i].step);
585
586     if (gmx_debug_at) {
587       fprintf(debug,"shell[%d] = %d\n",i,shell);
588       pr_rvec(debug,0,"fshell",f[shell],DIM,TRUE);
589       pr_rvec(debug,0,"xold",xcur[shell],DIM,TRUE);
590       pr_rvec(debug,0,"step",s[i].step,DIM,TRUE);
591       pr_rvec(debug,0,"xnew",xnew[shell],DIM,TRUE);
592     }
593   }
594 #ifdef PRINT_STEP
595   printf("step %.3e %.3e\n",step_min,step_max);
596 #endif
597 }
598
599 static void decrease_step_size(int nshell,t_shell s[])
600 {
601   int i;
602   
603   for(i=0; i<nshell; i++)
604     svmul(0.8,s[i].step,s[i].step);
605 }
606
607 static void print_epot(FILE *fp,gmx_large_int_t mdstep,int count,real epot,real df,
608                        int ndir,real sf_dir)
609 {
610   char buf[22];
611
612   fprintf(fp,"MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
613           gmx_step_str(mdstep,buf),count,epot,df);
614   if (ndir)
615     fprintf(fp,", dir. rmsF: %6.2e\n",sqrt(sf_dir/ndir));
616   else
617     fprintf(fp,"\n");
618 }
619
620
621 static real rms_force(t_commrec *cr,rvec f[],int ns,t_shell s[],
622                       int ndir,real *sf_dir,real *Epot)
623 {
624   int  i,shell,ntot;
625   double buf[4];
626
627   buf[0] = *sf_dir;
628   for(i=0; i<ns; i++) {
629     shell = s[i].shell;
630     buf[0]  += norm2(f[shell]);
631   }
632   ntot = ns;
633
634   if (PAR(cr)) {
635     buf[1] = ntot;
636     buf[2] = *sf_dir;
637     buf[3] = *Epot;
638     gmx_sumd(4,buf,cr);
639     ntot = (int)(buf[1] + 0.5);
640     *sf_dir = buf[2];
641     *Epot   = buf[3];
642   }
643   ntot += ndir;
644
645   return (ntot ? sqrt(buf[0]/ntot) : 0);
646 }
647
648 static void check_pbc(FILE *fp,rvec x[],int shell)
649 {
650   int m,now;
651   
652   now = shell-4;
653   for(m=0; (m<DIM); m++)
654     if (fabs(x[shell][m]-x[now][m]) > 0.3) {
655       pr_rvecs(fp,0,"SHELL-X",x+now,5);
656       break;
657     }
658 }
659
660 static void dump_shells(FILE *fp,rvec x[],rvec f[],real ftol,int ns,t_shell s[])
661 {
662   int  i,shell;
663   real ft2,ff2;
664   
665   ft2 = sqr(ftol);
666   
667   for(i=0; (i<ns); i++) {
668     shell = s[i].shell;
669     ff2   = iprod(f[shell],f[shell]);
670     if (ff2 > ft2)
671       fprintf(fp,"SHELL %5d, force %10.5f  %10.5f  %10.5f, |f| %10.5f\n",
672               shell,f[shell][XX],f[shell][YY],f[shell][ZZ],sqrt(ff2));
673     check_pbc(fp,x,shell);
674   }
675 }
676
677 static void init_adir(FILE *log,gmx_shellfc_t shfc,
678                       gmx_constr_t constr,t_idef *idef,t_inputrec *ir,
679                       t_commrec *cr,int dd_ac1,
680                       gmx_large_int_t step,t_mdatoms *md,int start,int end,
681                       rvec *x_old,rvec *x_init,rvec *x,
682                       rvec *f,rvec *acc_dir,matrix box,
683                       real lambda,real *dvdlambda,t_nrnb *nrnb)
684 {
685   rvec   *xnold,*xnew;
686   double w_dt;
687   int    gf,ga,gt;
688   real   dt,scale;
689   int    n,d; 
690   unsigned short *ptype;
691   rvec   p,dx;
692   
693   if (DOMAINDECOMP(cr))
694     n = dd_ac1;
695   else
696     n = end - start;
697   if (n > shfc->adir_nalloc) {
698     shfc->adir_nalloc = over_alloc_dd(n);
699     srenew(shfc->adir_xnold,shfc->adir_nalloc);
700     srenew(shfc->adir_xnew ,shfc->adir_nalloc);
701   }
702   xnold = shfc->adir_xnold;
703   xnew  = shfc->adir_xnew;
704     
705   ptype = md->ptype;
706
707   dt = ir->delta_t;
708
709   /* Does NOT work with freeze or acceleration groups (yet) */
710   for (n=start; n<end; n++) {  
711     w_dt = md->invmass[n]*dt;
712     
713     for (d=0; d<DIM; d++) {
714       if ((ptype[n] != eptVSite) && (ptype[n] != eptShell)) {
715         xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
716         xnew[n-start][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
717       } else {
718         xnold[n-start][d] = x[n][d];
719         xnew[n-start][d] = x[n][d];
720       }
721     }
722   }
723   constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
724             x,xnold-start,NULL,box,
725             lambda,dvdlambda,NULL,NULL,nrnb,econqCoord,FALSE,0,0);
726   constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
727             x,xnew-start,NULL,box,
728             lambda,dvdlambda,NULL,NULL,nrnb,econqCoord,FALSE,0,0);
729
730   /* Set xnew to minus the acceleration */
731   for (n=start; n<end; n++) {
732     for(d=0; d<DIM; d++)
733       xnew[n-start][d] =
734         -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
735         - f[n][d]*md->invmass[n];
736     clear_rvec(acc_dir[n]);
737   }
738
739   /* Project the acceleration on the old bond directions */
740   constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
741             x_old,xnew-start,acc_dir,box,
742             lambda,dvdlambda,NULL,NULL,nrnb,econqDeriv_FlexCon,FALSE,0,0); 
743 }
744
745 int relax_shell_flexcon(FILE *fplog,t_commrec *cr,gmx_bool bVerbose,
746                         gmx_large_int_t mdstep,t_inputrec *inputrec,
747                         gmx_bool bDoNS,int force_flags,
748                         gmx_bool bStopCM,
749                         gmx_localtop_t *top,
750                         gmx_mtop_t* mtop,
751                         gmx_constr_t constr,
752                         gmx_enerdata_t *enerd,t_fcdata *fcd,
753                         t_state *state,rvec f[],
754                         tensor force_vir,
755                         t_mdatoms *md,
756                         t_nrnb *nrnb,gmx_wallcycle_t wcycle,
757                         t_graph *graph,
758                         gmx_groups_t *groups,
759                         struct gmx_shellfc *shfc,
760                         t_forcerec *fr,
761                         gmx_bool bBornRadii,
762                         double t,rvec mu_tot,
763                         int natoms,gmx_bool *bConverged,
764                         gmx_vsite_t *vsite,
765                         FILE *fp_field)
766 {
767   int    nshell;
768   t_shell *shell;
769   t_idef *idef;
770   rvec   *pos[2],*force[2],*acc_dir=NULL,*x_old=NULL;
771   real   Epot[2],df[2];
772   rvec   dx;
773   real   sf_dir,invdt;
774   real   ftol,xiH,xiS,dum=0;
775   char   sbuf[22];
776   gmx_bool   bCont,bInit;
777   int    nat,dd_ac0,dd_ac1=0,i;
778   int    start=md->start,homenr=md->homenr,end=start+homenr,cg0,cg1;
779   int    nflexcon,g,number_steps,d,Min=0,count=0;
780 #define  Try (1-Min)             /* At start Try = 1 */
781
782   bCont        = (mdstep == inputrec->init_step) && inputrec->bContinuation;
783   bInit        = (mdstep == inputrec->init_step) || shfc->bForceInit;
784   ftol         = inputrec->em_tol;
785   number_steps = inputrec->niter;
786   nshell       = shfc->nshell;
787   shell        = shfc->shell;
788   nflexcon     = shfc->nflexcon;
789
790   idef = &top->idef;
791
792   if (DOMAINDECOMP(cr)) {
793     nat = dd_natoms_vsite(cr->dd);
794     if (nflexcon > 0) {
795       dd_get_constraint_range(cr->dd,&dd_ac0,&dd_ac1);
796       nat = max(nat,dd_ac1);
797     }
798   } else {
799     nat = state->natoms;
800   }
801
802   if (nat > shfc->x_nalloc) {
803     /* Allocate local arrays */
804     shfc->x_nalloc = over_alloc_dd(nat);
805     for(i=0; (i<2); i++) {
806       srenew(shfc->x[i],shfc->x_nalloc);
807       srenew(shfc->f[i],shfc->x_nalloc);
808     }
809   }
810   for(i=0; (i<2); i++) {
811     pos[i]   = shfc->x[i];
812     force[i] = shfc->f[i];
813   }
814      
815   /* With particle decomposition this code only works
816    * when all particles involved with each shell are in the same cg.
817    */
818
819   if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr)) {
820     /* This is the only time where the coordinates are used
821      * before do_force is called, which normally puts all
822      * charge groups in the box.
823      */
824     if (PARTDECOMP(cr)) {
825       pd_cg_range(cr,&cg0,&cg1);
826     } else {
827       cg0 = 0;
828       cg1 = top->cgs.nr;
829     }
830     put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,state->box,
831                              &(top->cgs),state->x,fr->cg_cm);
832     if (graph)
833       mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
834   }
835
836   /* After this all coordinate arrays will contain whole molecules */
837   if (graph)
838     shift_self(graph,state->box,state->x);
839
840   if (nflexcon) {
841     if (nat > shfc->flex_nalloc) {
842       shfc->flex_nalloc = over_alloc_dd(nat);
843       srenew(shfc->acc_dir,shfc->flex_nalloc);
844       srenew(shfc->x_old,shfc->flex_nalloc);
845     }
846     acc_dir = shfc->acc_dir;
847     x_old   = shfc->x_old;
848     for(i=0; i<homenr; i++) {
849       for(d=0; d<DIM; d++)
850         shfc->x_old[i][d] =
851           state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
852     }
853   }
854
855   /* Do a prediction of the shell positions */
856   if (shfc->bPredict && !bCont) {
857     predict_shells(fplog,state->x,state->v,inputrec->delta_t,nshell,shell,
858                    md->massT,NULL,bInit);
859   }
860
861   /* do_force expected the charge groups to be in the box */
862   if (graph)
863     unshift_self(graph,state->box,state->x);
864
865   /* Calculate the forces first time around */
866   if (gmx_debug_at) {
867     pr_rvecs(debug,0,"x b4 do_force",state->x + start,homenr);
868   }
869   do_force(fplog,cr,inputrec,mdstep,nrnb,wcycle,top,mtop,groups,
870            state->box,state->x,&state->hist,
871            force[Min],force_vir,md,enerd,fcd,
872            state->lambda,graph,
873            fr,vsite,mu_tot,t,fp_field,NULL,bBornRadii,
874            (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
875
876   sf_dir = 0;
877   if (nflexcon) {
878     init_adir(fplog,shfc,
879               constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
880               shfc->x_old-start,state->x,state->x,force[Min],
881               shfc->acc_dir-start,state->box,state->lambda,&dum,nrnb);
882
883     for(i=start; i<end; i++)
884       sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
885   }
886
887   Epot[Min] = enerd->term[F_EPOT];
888
889   df[Min]=rms_force(cr,shfc->f[Min],nshell,shell,nflexcon,&sf_dir,&Epot[Min]);
890   df[Try]=0;
891   if (debug) {
892     fprintf(debug,"df = %g  %g\n",df[Min],df[Try]);
893   }
894
895   if (gmx_debug_at) {
896     pr_rvecs(debug,0,"force0",force[Min],md->nr);
897   }
898
899   if (nshell+nflexcon > 0) {
900     /* Copy x to pos[Min] & pos[Try]: during minimization only the
901      * shell positions are updated, therefore the other particles must
902      * be set here.
903      */
904     memcpy(pos[Min],state->x,nat*sizeof(state->x[0]));
905     memcpy(pos[Try],state->x,nat*sizeof(state->x[0]));
906   }
907   
908   if (bVerbose && MASTER(cr))
909     print_epot(stdout,mdstep,0,Epot[Min],df[Min],nflexcon,sf_dir);
910
911   if (debug) {
912     fprintf(debug,"%17s: %14.10e\n",
913             interaction_function[F_EKIN].longname,enerd->term[F_EKIN]);
914     fprintf(debug,"%17s: %14.10e\n",
915             interaction_function[F_EPOT].longname,enerd->term[F_EPOT]);
916     fprintf(debug,"%17s: %14.10e\n",
917             interaction_function[F_ETOT].longname,enerd->term[F_ETOT]);
918     fprintf(debug,"SHELLSTEP %s\n",gmx_step_str(mdstep,sbuf));
919   }
920   
921   /* First check whether we should do shells, or whether the force is 
922    * low enough even without minimization.
923    */
924   *bConverged = (df[Min] < ftol);
925   
926   for(count=1; (!(*bConverged) && (count < number_steps)); count++) {
927     if (vsite)
928       construct_vsites(fplog,vsite,pos[Min],nrnb,inputrec->delta_t,state->v,
929                        idef->iparams,idef->il,
930                        fr->ePBC,fr->bMolPBC,graph,cr,state->box);
931      
932     if (nflexcon) {
933       init_adir(fplog,shfc,
934                 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
935                 x_old-start,state->x,pos[Min],force[Min],acc_dir-start,
936                 state->box,state->lambda,&dum,nrnb);
937       
938       directional_sd(fplog,pos[Min],pos[Try],acc_dir-start,start,end,
939                      fr->fc_stepsize);
940     }
941     
942     /* New positions, Steepest descent */
943     shell_pos_sd(fplog,pos[Min],pos[Try],force[Min],nshell,shell,count); 
944
945     /* do_force expected the charge groups to be in the box */
946     if (graph)
947       unshift_self(graph,state->box,pos[Try]);
948
949     if (gmx_debug_at) {
950       pr_rvecs(debug,0,"RELAX: pos[Min]  ",pos[Min] + start,homenr);
951       pr_rvecs(debug,0,"RELAX: pos[Try]  ",pos[Try] + start,homenr);
952     }
953     /* Try the new positions */
954     do_force(fplog,cr,inputrec,1,nrnb,wcycle,
955              top,mtop,groups,state->box,pos[Try],&state->hist,
956              force[Try],force_vir,
957              md,enerd,fcd,state->lambda,graph,
958              fr,vsite,mu_tot,t,fp_field,NULL,bBornRadii,
959              force_flags);
960     
961     if (gmx_debug_at) {
962       pr_rvecs(debug,0,"RELAX: force[Min]",force[Min] + start,homenr);
963       pr_rvecs(debug,0,"RELAX: force[Try]",force[Try] + start,homenr);
964     }
965     sf_dir = 0;
966     if (nflexcon) {
967       init_adir(fplog,shfc,
968                 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
969                 x_old-start,state->x,pos[Try],force[Try],acc_dir-start,
970                 state->box,state->lambda,&dum,nrnb);
971
972       for(i=start; i<end; i++)
973         sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
974     }
975
976     Epot[Try] = enerd->term[F_EPOT]; 
977     
978     df[Try]=rms_force(cr,force[Try],nshell,shell,nflexcon,&sf_dir,&Epot[Try]);
979
980     if (debug)
981       fprintf(debug,"df = %g  %g\n",df[Min],df[Try]);
982
983     if (debug) {
984       if (gmx_debug_at)
985         pr_rvecs(debug,0,"F na do_force",force[Try] + start,homenr);
986       if (gmx_debug_at) {
987         fprintf(debug,"SHELL ITER %d\n",count);
988         dump_shells(debug,pos[Try],force[Try],ftol,nshell,shell);
989       }
990     }
991
992     if (bVerbose && MASTER(cr))
993       print_epot(stdout,mdstep,count,Epot[Try],df[Try],nflexcon,sf_dir);
994       
995     *bConverged = (df[Try] < ftol);
996     
997     if ((df[Try] < df[Min])) {
998       if (debug)
999         fprintf(debug,"Swapping Min and Try\n");
1000       if (nflexcon) {
1001         /* Correct the velocities for the flexible constraints */
1002         invdt = 1/inputrec->delta_t;
1003         for(i=start; i<end; i++) {
1004           for(d=0; d<DIM; d++)
1005             state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1006         }
1007       }
1008       Min  = Try;
1009     } else {
1010       decrease_step_size(nshell,shell);
1011     }
1012   }
1013   if (MASTER(cr) && !(*bConverged)) {
1014     /* Note that the energies and virial are incorrect when not converged */
1015     if (fplog)
1016       fprintf(fplog,
1017               "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1018               gmx_step_str(mdstep,sbuf),number_steps,df[Min]);
1019     fprintf(stderr,
1020             "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1021             gmx_step_str(mdstep,sbuf),number_steps,df[Min]);
1022   }
1023
1024   /* Copy back the coordinates and the forces */
1025   memcpy(state->x,pos[Min],nat*sizeof(state->x[0]));
1026   memcpy(f,force[Min],nat*sizeof(f[0]));
1027
1028   return count; 
1029 }
1030