Merge branch 'release-4-6'
[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 "force.h"
45 #include "mdrun.h"
46 #include "partdec.h"
47 #include "mdatoms.h"
48 #include "vsite.h"
49 #include "network.h"
50 #include "names.h"
51 #include "constr.h"
52 #include "domdec.h"
53 #include "partdec.h"
54 #include "physics.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    bRequireInit;   /* Require 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     {
105         fprintf(fplog, "%5d  %8.3f  %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
106         if (s[i].nnucl == 2)
107         {
108             fprintf(fplog, "  %5d\n", s[i].nucl2);
109         }
110         else if (s[i].nnucl == 3)
111         {
112             fprintf(fplog, "  %5d  %5d\n", s[i].nucl2, s[i].nucl3);
113         }
114         else
115         {
116             fprintf(fplog, "\n");
117         }
118     }
119 }
120
121 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
122                            int ns, t_shell s[],
123                            real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
124 {
125     int                   i, m, s1, n1, n2, n3;
126     real                  dt_1, dt_2, dt_3, fudge, tm, m1, m2, m3;
127     rvec                 *ptr;
128     gmx_mtop_atomlookup_t alook = NULL;
129     t_atom               *atom;
130
131     if (mass == NULL)
132     {
133         alook = gmx_mtop_atomlookup_init(mtop);
134     }
135
136     /* We introduce a fudge factor for performance reasons: with this choice
137      * the initial force on the shells is about a factor of two lower than
138      * without
139      */
140     fudge = 1.0;
141
142     if (bInit)
143     {
144         if (fplog)
145         {
146             fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
147         }
148         ptr  = x;
149         dt_1 = 1;
150     }
151     else
152     {
153         ptr  = v;
154         dt_1 = fudge*dt;
155     }
156
157     for (i = 0; (i < ns); i++)
158     {
159         s1 = s[i].shell;
160         if (bInit)
161         {
162             clear_rvec(x[s1]);
163         }
164         switch (s[i].nnucl)
165         {
166             case 1:
167                 n1 = s[i].nucl1;
168                 for (m = 0; (m < DIM); m++)
169                 {
170                     x[s1][m] += ptr[n1][m]*dt_1;
171                 }
172                 break;
173             case 2:
174                 n1 = s[i].nucl1;
175                 n2 = s[i].nucl2;
176                 if (mass)
177                 {
178                     m1 = mass[n1];
179                     m2 = mass[n2];
180                 }
181                 else
182                 {
183                     /* Not the correct masses with FE, but it is just a prediction... */
184                     m1 = atom[n1].m;
185                     m2 = atom[n2].m;
186                 }
187                 tm = dt_1/(m1+m2);
188                 for (m = 0; (m < DIM); m++)
189                 {
190                     x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
191                 }
192                 break;
193             case 3:
194                 n1 = s[i].nucl1;
195                 n2 = s[i].nucl2;
196                 n3 = s[i].nucl3;
197                 if (mass)
198                 {
199                     m1 = mass[n1];
200                     m2 = mass[n2];
201                     m3 = mass[n3];
202                 }
203                 else
204                 {
205                     /* Not the correct masses with FE, but it is just a prediction... */
206                     gmx_mtop_atomnr_to_atom(alook, n1, &atom);
207                     m1 = atom->m;
208                     gmx_mtop_atomnr_to_atom(alook, n2, &atom);
209                     m2 = atom->m;
210                     gmx_mtop_atomnr_to_atom(alook, n3, &atom);
211                     m3 = atom->m;
212                 }
213                 tm = dt_1/(m1+m2+m3);
214                 for (m = 0; (m < DIM); m++)
215                 {
216                     x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
217                 }
218                 break;
219             default:
220                 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
221         }
222     }
223
224     if (mass == NULL)
225     {
226         gmx_mtop_atomlookup_destroy(alook);
227     }
228 }
229
230 gmx_shellfc_t init_shell_flexcon(FILE *fplog,
231                                  gmx_mtop_t *mtop, int nflexcon,
232                                  rvec *x)
233 {
234     struct gmx_shellfc       *shfc;
235     t_shell                  *shell;
236     int                      *shell_index = NULL, *at2cg;
237     t_atom                   *atom;
238     int                       n[eptNR], ns, nshell, nsi;
239     int                       i, j, nmol, type, mb, mt, a_offset, cg, mol, ftype, nra;
240     real                      qS, alpha;
241     int                       aS, aN = 0; /* Shell and nucleus */
242     int                       bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
243 #define NBT asize(bondtypes)
244     t_iatom                  *ia;
245     gmx_mtop_atomloop_block_t aloopb;
246     gmx_mtop_atomloop_all_t   aloop;
247     gmx_ffparams_t           *ffparams;
248     gmx_molblock_t           *molb;
249     gmx_moltype_t            *molt;
250     t_block                  *cgs;
251
252     /* Count number of shells, and find their indices */
253     for (i = 0; (i < eptNR); i++)
254     {
255         n[i] = 0;
256     }
257
258     aloopb = gmx_mtop_atomloop_block_init(mtop);
259     while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
260     {
261         n[atom->ptype] += nmol;
262     }
263
264     if (fplog)
265     {
266         /* Print the number of each particle type */
267         for (i = 0; (i < eptNR); i++)
268         {
269             if (n[i] != 0)
270             {
271                 fprintf(fplog, "There are: %d %ss\n", n[i], ptype_str[i]);
272             }
273         }
274     }
275
276     nshell = n[eptShell];
277
278     if (nshell == 0 && nflexcon == 0)
279     {
280         return NULL;
281     }
282
283     snew(shfc, 1);
284     shfc->nflexcon = nflexcon;
285
286     if (nshell == 0)
287     {
288         return shfc;
289     }
290
291     /* We have shells: fill the shell data structure */
292
293     /* Global system sized array, this should be avoided */
294     snew(shell_index, mtop->natoms);
295
296     aloop  = gmx_mtop_atomloop_all_init(mtop);
297     nshell = 0;
298     while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
299     {
300         if (atom->ptype == eptShell)
301         {
302             shell_index[i] = nshell++;
303         }
304     }
305
306     snew(shell, nshell);
307
308     /* Initiate the shell structures */
309     for (i = 0; (i < nshell); i++)
310     {
311         shell[i].shell = NO_ATID;
312         shell[i].nnucl = 0;
313         shell[i].nucl1 = NO_ATID;
314         shell[i].nucl2 = NO_ATID;
315         shell[i].nucl3 = NO_ATID;
316         /* shell[i].bInterCG=FALSE; */
317         shell[i].k_1   = 0;
318         shell[i].k     = 0;
319     }
320
321     ffparams = &mtop->ffparams;
322
323     /* Now fill the structures */
324     shfc->bInterCG = FALSE;
325     ns             = 0;
326     a_offset       = 0;
327     for (mb = 0; mb < mtop->nmolblock; mb++)
328     {
329         molb = &mtop->molblock[mb];
330         molt = &mtop->moltype[molb->type];
331
332         cgs = &molt->cgs;
333         snew(at2cg, molt->atoms.nr);
334         for (cg = 0; cg < cgs->nr; cg++)
335         {
336             for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
337             {
338                 at2cg[i] = cg;
339             }
340         }
341
342         atom = molt->atoms.atom;
343         for (mol = 0; mol < molb->nmol; mol++)
344         {
345             for (j = 0; (j < NBT); j++)
346             {
347                 ia = molt->ilist[bondtypes[j]].iatoms;
348                 for (i = 0; (i < molt->ilist[bondtypes[j]].nr); )
349                 {
350                     type  = ia[0];
351                     ftype = ffparams->functype[type];
352                     nra   = interaction_function[ftype].nratoms;
353
354                     /* Check whether we have a bond with a shell */
355                     aS = NO_ATID;
356
357                     switch (bondtypes[j])
358                     {
359                         case F_BONDS:
360                         case F_HARMONIC:
361                         case F_CUBICBONDS:
362                         case F_POLARIZATION:
363                         case F_ANHARM_POL:
364                             if (atom[ia[1]].ptype == eptShell)
365                             {
366                                 aS = ia[1];
367                                 aN = ia[2];
368                             }
369                             else if (atom[ia[2]].ptype == eptShell)
370                             {
371                                 aS = ia[2];
372                                 aN = ia[1];
373                             }
374                             break;
375                         case F_WATER_POL:
376                             aN    = ia[4]; /* Dummy */
377                             aS    = ia[5]; /* Shell */
378                             break;
379                         default:
380                             gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
381                     }
382
383                     if (aS != NO_ATID)
384                     {
385                         qS = atom[aS].q;
386
387                         /* Check whether one of the particles is a shell... */
388                         nsi = shell_index[a_offset+aS];
389                         if ((nsi < 0) || (nsi >= nshell))
390                         {
391                             gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
392                                       nsi, nshell, aS);
393                         }
394                         if (shell[nsi].shell == NO_ATID)
395                         {
396                             shell[nsi].shell = a_offset + aS;
397                             ns++;
398                         }
399                         else if (shell[nsi].shell != a_offset+aS)
400                         {
401                             gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
402                         }
403
404                         if      (shell[nsi].nucl1 == NO_ATID)
405                         {
406                             shell[nsi].nucl1 = a_offset + aN;
407                         }
408                         else if (shell[nsi].nucl2 == NO_ATID)
409                         {
410                             shell[nsi].nucl2 = a_offset + aN;
411                         }
412                         else if (shell[nsi].nucl3 == NO_ATID)
413                         {
414                             shell[nsi].nucl3 = a_offset + aN;
415                         }
416                         else
417                         {
418                             if (fplog)
419                             {
420                                 pr_shell(fplog, ns, shell);
421                             }
422                             gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
423                         }
424                         if (at2cg[aS] != at2cg[aN])
425                         {
426                             /* shell[nsi].bInterCG = TRUE; */
427                             shfc->bInterCG = TRUE;
428                         }
429
430                         switch (bondtypes[j])
431                         {
432                             case F_BONDS:
433                             case F_HARMONIC:
434                                 shell[nsi].k    += ffparams->iparams[type].harmonic.krA;
435                                 break;
436                             case F_CUBICBONDS:
437                                 shell[nsi].k    += ffparams->iparams[type].cubic.kb;
438                                 break;
439                             case F_POLARIZATION:
440                             case F_ANHARM_POL:
441                                 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
442                                 {
443                                     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);
444                                 }
445                                 shell[nsi].k    += sqr(qS)*ONE_4PI_EPS0/
446                                     ffparams->iparams[type].polarize.alpha;
447                                 break;
448                             case F_WATER_POL:
449                                 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
450                                 {
451                                     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);
452                                 }
453                                 alpha          = (ffparams->iparams[type].wpol.al_x+
454                                                   ffparams->iparams[type].wpol.al_y+
455                                                   ffparams->iparams[type].wpol.al_z)/3.0;
456                                 shell[nsi].k  += sqr(qS)*ONE_4PI_EPS0/alpha;
457                                 break;
458                             default:
459                                 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
460                         }
461                         shell[nsi].nnucl++;
462                     }
463                     ia += nra+1;
464                     i  += nra+1;
465                 }
466             }
467             a_offset += molt->atoms.nr;
468         }
469         /* Done with this molecule type */
470         sfree(at2cg);
471     }
472
473     /* Verify whether it's all correct */
474     if (ns != nshell)
475     {
476         gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
477     }
478
479     for (i = 0; (i < ns); i++)
480     {
481         shell[i].k_1 = 1.0/shell[i].k;
482     }
483
484     if (debug)
485     {
486         pr_shell(debug, ns, shell);
487     }
488
489
490     shfc->nshell_gl      = ns;
491     shfc->shell_gl       = shell;
492     shfc->shell_index_gl = shell_index;
493
494     shfc->bPredict     = (getenv("GMX_NOPREDICT") == NULL);
495     shfc->bRequireInit = FALSE;
496     if (!shfc->bPredict)
497     {
498         if (fplog)
499         {
500             fprintf(fplog, "\nWill never predict shell positions\n");
501         }
502     }
503     else
504     {
505         shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != NULL);
506         if (shfc->bRequireInit && fplog)
507         {
508             fprintf(fplog, "\nWill always initiate shell positions\n");
509         }
510     }
511
512     if (shfc->bPredict)
513     {
514         if (x)
515         {
516             predict_shells(fplog, x, NULL, 0, shfc->nshell_gl, shfc->shell_gl,
517                            NULL, mtop, TRUE);
518         }
519
520         if (shfc->bInterCG)
521         {
522             if (fplog)
523             {
524                 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");
525             }
526             shfc->bPredict = FALSE;
527         }
528     }
529
530     return shfc;
531 }
532
533 void make_local_shells(t_commrec *cr, t_mdatoms *md,
534                        struct gmx_shellfc *shfc)
535 {
536     t_shell      *shell;
537     int           a0, a1, *ind, nshell, i;
538     gmx_domdec_t *dd = NULL;
539
540     if (PAR(cr))
541     {
542         if (DOMAINDECOMP(cr))
543         {
544             dd = cr->dd;
545             a0 = 0;
546             a1 = dd->nat_home;
547         }
548         else
549         {
550             pd_at_range(cr, &a0, &a1);
551         }
552     }
553     else
554     {
555         /* Single node: we need all shells, just copy the pointer */
556         shfc->nshell = shfc->nshell_gl;
557         shfc->shell  = shfc->shell_gl;
558
559         return;
560     }
561
562     ind = shfc->shell_index_gl;
563
564     nshell = 0;
565     shell  = shfc->shell;
566     for (i = a0; i < a1; i++)
567     {
568         if (md->ptype[i] == eptShell)
569         {
570             if (nshell+1 > shfc->shell_nalloc)
571             {
572                 shfc->shell_nalloc = over_alloc_dd(nshell+1);
573                 srenew(shell, shfc->shell_nalloc);
574             }
575             if (dd)
576             {
577                 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
578             }
579             else
580             {
581                 shell[nshell] = shfc->shell_gl[ind[i]];
582             }
583             /* With inter-cg shells we can no do shell prediction,
584              * so we do not need the nuclei numbers.
585              */
586             if (!shfc->bInterCG)
587             {
588                 shell[nshell].nucl1   = i + shell[nshell].nucl1 - shell[nshell].shell;
589                 if (shell[nshell].nnucl > 1)
590                 {
591                     shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
592                 }
593                 if (shell[nshell].nnucl > 2)
594                 {
595                     shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
596                 }
597             }
598             shell[nshell].shell = i;
599             nshell++;
600         }
601     }
602
603     shfc->nshell = nshell;
604     shfc->shell  = shell;
605 }
606
607 static void do_1pos(rvec xnew, rvec xold, rvec f, real step)
608 {
609     real xo, yo, zo;
610     real dx, dy, dz;
611
612     xo = xold[XX];
613     yo = xold[YY];
614     zo = xold[ZZ];
615
616     dx = f[XX]*step;
617     dy = f[YY]*step;
618     dz = f[ZZ]*step;
619
620     xnew[XX] = xo+dx;
621     xnew[YY] = yo+dy;
622     xnew[ZZ] = zo+dz;
623 }
624
625 static void do_1pos3(rvec xnew, rvec xold, rvec f, rvec step)
626 {
627     real xo, yo, zo;
628     real dx, dy, dz;
629
630     xo = xold[XX];
631     yo = xold[YY];
632     zo = xold[ZZ];
633
634     dx = f[XX]*step[XX];
635     dy = f[YY]*step[YY];
636     dz = f[ZZ]*step[ZZ];
637
638     xnew[XX] = xo+dx;
639     xnew[YY] = yo+dy;
640     xnew[ZZ] = zo+dz;
641 }
642
643 static void directional_sd(rvec xold[], rvec xnew[], rvec acc_dir[],
644                            int start, int homenr, real step)
645 {
646     int  i;
647
648     for (i = start; i < homenr; i++)
649     {
650         do_1pos(xnew[i], xold[i], acc_dir[i], step);
651     }
652 }
653
654 static void shell_pos_sd(rvec xcur[], rvec xnew[], rvec f[],
655                          int ns, t_shell s[], int count)
656 {
657     const real step_scale_min       = 0.8,
658                step_scale_increment = 0.2,
659                step_scale_max       = 1.2,
660                step_scale_multiple  = (step_scale_max - step_scale_min) / step_scale_increment;
661     int  i, shell, d;
662     real dx, df, k_est;
663 #ifdef PRINT_STEP
664     real step_min, step_max;
665
666     step_min = 1e30;
667     step_max = 0;
668 #endif
669     for (i = 0; (i < ns); i++)
670     {
671         shell = s[i].shell;
672         if (count == 1)
673         {
674             for (d = 0; d < DIM; d++)
675             {
676                 s[i].step[d] = s[i].k_1;
677 #ifdef PRINT_STEP
678                 step_min = min(step_min, s[i].step[d]);
679                 step_max = max(step_max, s[i].step[d]);
680 #endif
681             }
682         }
683         else
684         {
685             for (d = 0; d < DIM; d++)
686             {
687                 dx = xcur[shell][d] - s[i].xold[d];
688                 df =    f[shell][d] - s[i].fold[d];
689                 /* -dx/df gets used to generate an interpolated value, but would
690                  * cause a NaN if df were binary-equal to zero. Values close to
691                  * zero won't cause problems (because of the min() and max()), so
692                  * just testing for binary inequality is OK. */
693                 if (0.0 != df)
694                 {
695                     k_est = -dx/df;
696                     /* Scale the step size by a factor interpolated from
697                      * step_scale_min to step_scale_max, as k_est goes from 0 to
698                      * step_scale_multiple * s[i].step[d] */
699                     s[i].step[d] =
700                         step_scale_min * s[i].step[d] +
701                         step_scale_increment * min(step_scale_multiple * s[i].step[d], max(k_est, 0));
702                 }
703                 else
704                 {
705                     /* Here 0 == df */
706                     if (gmx_numzero(dx)) /* 0 == dx */
707                     {
708                         /* Likely this will never happen, but if it does just
709                          * don't scale the step. */
710                     }
711                     else /* 0 != dx */
712                     {
713                         s[i].step[d] *= step_scale_max;
714                     }
715                 }
716 #ifdef PRINT_STEP
717                 step_min = min(step_min, s[i].step[d]);
718                 step_max = max(step_max, s[i].step[d]);
719 #endif
720             }
721         }
722         copy_rvec(xcur[shell], s[i].xold);
723         copy_rvec(f[shell],   s[i].fold);
724
725         do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
726
727         if (gmx_debug_at)
728         {
729             fprintf(debug, "shell[%d] = %d\n", i, shell);
730             pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
731             pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
732             pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
733             pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
734         }
735     }
736 #ifdef PRINT_STEP
737     printf("step %.3e %.3e\n", step_min, step_max);
738 #endif
739 }
740
741 static void decrease_step_size(int nshell, t_shell s[])
742 {
743     int i;
744
745     for (i = 0; i < nshell; i++)
746     {
747         svmul(0.8, s[i].step, s[i].step);
748     }
749 }
750
751 static void print_epot(FILE *fp, gmx_large_int_t mdstep, int count, real epot, real df,
752                        int ndir, real sf_dir)
753 {
754     char buf[22];
755
756     fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
757             gmx_step_str(mdstep, buf), count, epot, df);
758     if (ndir)
759     {
760         fprintf(fp, ", dir. rmsF: %6.2e\n", sqrt(sf_dir/ndir));
761     }
762     else
763     {
764         fprintf(fp, "\n");
765     }
766 }
767
768
769 static real rms_force(t_commrec *cr, rvec f[], int ns, t_shell s[],
770                       int ndir, real *sf_dir, real *Epot)
771 {
772     int    i, shell, ntot;
773     double buf[4];
774
775     buf[0] = *sf_dir;
776     for (i = 0; i < ns; i++)
777     {
778         shell    = s[i].shell;
779         buf[0]  += norm2(f[shell]);
780     }
781     ntot = ns;
782
783     if (PAR(cr))
784     {
785         buf[1] = ntot;
786         buf[2] = *sf_dir;
787         buf[3] = *Epot;
788         gmx_sumd(4, buf, cr);
789         ntot    = (int)(buf[1] + 0.5);
790         *sf_dir = buf[2];
791         *Epot   = buf[3];
792     }
793     ntot += ndir;
794
795     return (ntot ? sqrt(buf[0]/ntot) : 0);
796 }
797
798 static void check_pbc(FILE *fp, rvec x[], int shell)
799 {
800     int m, now;
801
802     now = shell-4;
803     for (m = 0; (m < DIM); m++)
804     {
805         if (fabs(x[shell][m]-x[now][m]) > 0.3)
806         {
807             pr_rvecs(fp, 0, "SHELL-X", x+now, 5);
808             break;
809         }
810     }
811 }
812
813 static void dump_shells(FILE *fp, rvec x[], rvec f[], real ftol, int ns, t_shell s[])
814 {
815     int  i, shell;
816     real ft2, ff2;
817
818     ft2 = sqr(ftol);
819
820     for (i = 0; (i < ns); i++)
821     {
822         shell = s[i].shell;
823         ff2   = iprod(f[shell], f[shell]);
824         if (ff2 > ft2)
825         {
826             fprintf(fp, "SHELL %5d, force %10.5f  %10.5f  %10.5f, |f| %10.5f\n",
827                     shell, f[shell][XX], f[shell][YY], f[shell][ZZ], sqrt(ff2));
828         }
829         check_pbc(fp, x, shell);
830     }
831 }
832
833 static void init_adir(FILE *log, gmx_shellfc_t shfc,
834                       gmx_constr_t constr, t_idef *idef, t_inputrec *ir,
835                       t_commrec *cr, int dd_ac1,
836                       gmx_large_int_t step, t_mdatoms *md, int start, int end,
837                       rvec *x_old, rvec *x_init, rvec *x,
838                       rvec *f, rvec *acc_dir,
839                       gmx_bool bMolPBC, matrix box,
840                       real *lambda, real *dvdlambda, t_nrnb *nrnb)
841 {
842     rvec           *xnold, *xnew;
843     double          w_dt;
844     int             gf, ga, gt;
845     real            dt, scale;
846     int             n, d;
847     unsigned short *ptype;
848     rvec            p, dx;
849
850     if (DOMAINDECOMP(cr))
851     {
852         n = dd_ac1;
853     }
854     else
855     {
856         n = end - start;
857     }
858     if (n > shfc->adir_nalloc)
859     {
860         shfc->adir_nalloc = over_alloc_dd(n);
861         srenew(shfc->adir_xnold, shfc->adir_nalloc);
862         srenew(shfc->adir_xnew, shfc->adir_nalloc);
863     }
864     xnold = shfc->adir_xnold;
865     xnew  = shfc->adir_xnew;
866
867     ptype = md->ptype;
868
869     dt = ir->delta_t;
870
871     /* Does NOT work with freeze or acceleration groups (yet) */
872     for (n = start; n < end; n++)
873     {
874         w_dt = md->invmass[n]*dt;
875
876         for (d = 0; d < DIM; d++)
877         {
878             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
879             {
880                 xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
881                 xnew[n-start][d]  = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
882             }
883             else
884             {
885                 xnold[n-start][d] = x[n][d];
886                 xnew[n-start][d]  = x[n][d];
887             }
888         }
889     }
890     constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
891               x, xnold-start, NULL, bMolPBC, box,
892               lambda[efptBONDED], &(dvdlambda[efptBONDED]),
893               NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
894     constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
895               x, xnew-start, NULL, bMolPBC, box,
896               lambda[efptBONDED], &(dvdlambda[efptBONDED]),
897               NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
898
899     for (n = start; n < end; n++)
900     {
901         for (d = 0; d < DIM; d++)
902         {
903             xnew[n-start][d] =
904                 -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
905                 - f[n][d]*md->invmass[n];
906         }
907         clear_rvec(acc_dir[n]);
908     }
909
910     /* Project the acceleration on the old bond directions */
911     constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
912               x_old, xnew-start, acc_dir, bMolPBC, box,
913               lambda[efptBONDED], &(dvdlambda[efptBONDED]),
914               NULL, NULL, nrnb, econqDeriv_FlexCon, FALSE, 0, 0);
915 }
916
917 int relax_shell_flexcon(FILE *fplog, t_commrec *cr, gmx_bool bVerbose,
918                         gmx_large_int_t mdstep, t_inputrec *inputrec,
919                         gmx_bool bDoNS, int force_flags,
920                         gmx_localtop_t *top,
921                         gmx_constr_t constr,
922                         gmx_enerdata_t *enerd, t_fcdata *fcd,
923                         t_state *state, rvec f[],
924                         tensor force_vir,
925                         t_mdatoms *md,
926                         t_nrnb *nrnb, gmx_wallcycle_t wcycle,
927                         t_graph *graph,
928                         gmx_groups_t *groups,
929                         struct gmx_shellfc *shfc,
930                         t_forcerec *fr,
931                         gmx_bool bBornRadii,
932                         double t, rvec mu_tot,
933                         gmx_bool *bConverged,
934                         gmx_vsite_t *vsite,
935                         FILE *fp_field)
936 {
937     int        nshell;
938     t_shell   *shell;
939     t_idef    *idef;
940     rvec      *pos[2], *force[2], *acc_dir = NULL, *x_old = NULL;
941     real       Epot[2], df[2];
942     rvec       dx;
943     real       sf_dir, invdt;
944     real       ftol, xiH, xiS, dum = 0;
945     char       sbuf[22];
946     gmx_bool   bCont, bInit;
947     int        nat, dd_ac0, dd_ac1 = 0, i;
948     int        start = md->start, homenr = md->homenr, end = start+homenr, cg0, cg1;
949     int        nflexcon, g, number_steps, d, Min = 0, count = 0;
950 #define  Try (1-Min)             /* At start Try = 1 */
951
952     bCont        = (mdstep == inputrec->init_step) && inputrec->bContinuation;
953     bInit        = (mdstep == inputrec->init_step) || shfc->bRequireInit;
954     ftol         = inputrec->em_tol;
955     number_steps = inputrec->niter;
956     nshell       = shfc->nshell;
957     shell        = shfc->shell;
958     nflexcon     = shfc->nflexcon;
959
960     idef = &top->idef;
961
962     if (DOMAINDECOMP(cr))
963     {
964         nat = dd_natoms_vsite(cr->dd);
965         if (nflexcon > 0)
966         {
967             dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
968             nat = max(nat, dd_ac1);
969         }
970     }
971     else
972     {
973         nat = state->natoms;
974     }
975
976     if (nat > shfc->x_nalloc)
977     {
978         /* Allocate local arrays */
979         shfc->x_nalloc = over_alloc_dd(nat);
980         for (i = 0; (i < 2); i++)
981         {
982             srenew(shfc->x[i], shfc->x_nalloc);
983             srenew(shfc->f[i], shfc->x_nalloc);
984         }
985     }
986     for (i = 0; (i < 2); i++)
987     {
988         pos[i]   = shfc->x[i];
989         force[i] = shfc->f[i];
990     }
991
992     /* With particle decomposition this code only works
993      * when all particles involved with each shell are in the same cg.
994      */
995
996     if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
997     {
998         /* This is the only time where the coordinates are used
999          * before do_force is called, which normally puts all
1000          * charge groups in the box.
1001          */
1002         if (PARTDECOMP(cr))
1003         {
1004             pd_cg_range(cr, &cg0, &cg1);
1005         }
1006         else
1007         {
1008             cg0 = 0;
1009             cg1 = top->cgs.nr;
1010         }
1011         put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1012                                  &(top->cgs), state->x, fr->cg_cm);
1013         if (graph)
1014         {
1015             mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
1016         }
1017     }
1018
1019     /* After this all coordinate arrays will contain whole molecules */
1020     if (graph)
1021     {
1022         shift_self(graph, state->box, state->x);
1023     }
1024
1025     if (nflexcon)
1026     {
1027         if (nat > shfc->flex_nalloc)
1028         {
1029             shfc->flex_nalloc = over_alloc_dd(nat);
1030             srenew(shfc->acc_dir, shfc->flex_nalloc);
1031             srenew(shfc->x_old, shfc->flex_nalloc);
1032         }
1033         acc_dir = shfc->acc_dir;
1034         x_old   = shfc->x_old;
1035         for (i = 0; i < homenr; i++)
1036         {
1037             for (d = 0; d < DIM; d++)
1038             {
1039                 shfc->x_old[i][d] =
1040                     state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
1041             }
1042         }
1043     }
1044
1045     /* Do a prediction of the shell positions */
1046     if (shfc->bPredict && !bCont)
1047     {
1048         predict_shells(fplog, state->x, state->v, inputrec->delta_t, nshell, shell,
1049                        md->massT, NULL, bInit);
1050     }
1051
1052     /* do_force expected the charge groups to be in the box */
1053     if (graph)
1054     {
1055         unshift_self(graph, state->box, state->x);
1056     }
1057
1058     /* Calculate the forces first time around */
1059     if (gmx_debug_at)
1060     {
1061         pr_rvecs(debug, 0, "x b4 do_force", state->x + start, homenr);
1062     }
1063     do_force(fplog, cr, inputrec, mdstep, nrnb, wcycle, top, groups,
1064              state->box, state->x, &state->hist,
1065              force[Min], force_vir, md, enerd, fcd,
1066              state->lambda, graph,
1067              fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1068              (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
1069
1070     sf_dir = 0;
1071     if (nflexcon)
1072     {
1073         init_adir(fplog, shfc,
1074                   constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1075                   shfc->x_old-start, state->x, state->x, force[Min],
1076                   shfc->acc_dir-start,
1077                   fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1078
1079         for (i = start; i < end; i++)
1080         {
1081             sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
1082         }
1083     }
1084
1085     Epot[Min] = enerd->term[F_EPOT];
1086
1087     df[Min] = rms_force(cr, shfc->f[Min], nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1088     df[Try] = 0;
1089     if (debug)
1090     {
1091         fprintf(debug, "df = %g  %g\n", df[Min], df[Try]);
1092     }
1093
1094     if (gmx_debug_at)
1095     {
1096         pr_rvecs(debug, 0, "force0", force[Min], md->nr);
1097     }
1098
1099     if (nshell+nflexcon > 0)
1100     {
1101         /* Copy x to pos[Min] & pos[Try]: during minimization only the
1102          * shell positions are updated, therefore the other particles must
1103          * be set here.
1104          */
1105         memcpy(pos[Min], state->x, nat*sizeof(state->x[0]));
1106         memcpy(pos[Try], state->x, nat*sizeof(state->x[0]));
1107     }
1108
1109     if (bVerbose && MASTER(cr))
1110     {
1111         print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1112     }
1113
1114     if (debug)
1115     {
1116         fprintf(debug, "%17s: %14.10e\n",
1117                 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1118         fprintf(debug, "%17s: %14.10e\n",
1119                 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1120         fprintf(debug, "%17s: %14.10e\n",
1121                 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1122         fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1123     }
1124
1125     /* First check whether we should do shells, or whether the force is
1126      * low enough even without minimization.
1127      */
1128     *bConverged = (df[Min] < ftol);
1129
1130     for (count = 1; (!(*bConverged) && (count < number_steps)); count++)
1131     {
1132         if (vsite)
1133         {
1134             construct_vsites(vsite, pos[Min], inputrec->delta_t, state->v,
1135                              idef->iparams, idef->il,
1136                              fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1137         }
1138
1139         if (nflexcon)
1140         {
1141             init_adir(fplog, shfc,
1142                       constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1143                       x_old-start, state->x, pos[Min], force[Min], acc_dir-start,
1144                       fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1145
1146             directional_sd(pos[Min], pos[Try], acc_dir-start, start, end,
1147                            fr->fc_stepsize);
1148         }
1149
1150         /* New positions, Steepest descent */
1151         shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1152
1153         /* do_force expected the charge groups to be in the box */
1154         if (graph)
1155         {
1156             unshift_self(graph, state->box, pos[Try]);
1157         }
1158
1159         if (gmx_debug_at)
1160         {
1161             pr_rvecs(debug, 0, "RELAX: pos[Min]  ", pos[Min] + start, homenr);
1162             pr_rvecs(debug, 0, "RELAX: pos[Try]  ", pos[Try] + start, homenr);
1163         }
1164         /* Try the new positions */
1165         do_force(fplog, cr, inputrec, 1, nrnb, wcycle,
1166                  top, groups, state->box, pos[Try], &state->hist,
1167                  force[Try], force_vir,
1168                  md, enerd, fcd, state->lambda, graph,
1169                  fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1170                  force_flags);
1171
1172         if (gmx_debug_at)
1173         {
1174             pr_rvecs(debug, 0, "RELAX: force[Min]", force[Min] + start, homenr);
1175             pr_rvecs(debug, 0, "RELAX: force[Try]", force[Try] + start, homenr);
1176         }
1177         sf_dir = 0;
1178         if (nflexcon)
1179         {
1180             init_adir(fplog, shfc,
1181                       constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1182                       x_old-start, state->x, pos[Try], force[Try], acc_dir-start,
1183                       fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1184
1185             for (i = start; i < end; i++)
1186             {
1187                 sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
1188             }
1189         }
1190
1191         Epot[Try] = enerd->term[F_EPOT];
1192
1193         df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1194
1195         if (debug)
1196         {
1197             fprintf(debug, "df = %g  %g\n", df[Min], df[Try]);
1198         }
1199
1200         if (debug)
1201         {
1202             if (gmx_debug_at)
1203             {
1204                 pr_rvecs(debug, 0, "F na do_force", force[Try] + start, homenr);
1205             }
1206             if (gmx_debug_at)
1207             {
1208                 fprintf(debug, "SHELL ITER %d\n", count);
1209                 dump_shells(debug, pos[Try], force[Try], ftol, nshell, shell);
1210             }
1211         }
1212
1213         if (bVerbose && MASTER(cr))
1214         {
1215             print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1216         }
1217
1218         *bConverged = (df[Try] < ftol);
1219
1220         if ((df[Try] < df[Min]))
1221         {
1222             if (debug)
1223             {
1224                 fprintf(debug, "Swapping Min and Try\n");
1225             }
1226             if (nflexcon)
1227             {
1228                 /* Correct the velocities for the flexible constraints */
1229                 invdt = 1/inputrec->delta_t;
1230                 for (i = start; i < end; i++)
1231                 {
1232                     for (d = 0; d < DIM; d++)
1233                     {
1234                         state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1235                     }
1236                 }
1237             }
1238             Min  = Try;
1239         }
1240         else
1241         {
1242             decrease_step_size(nshell, shell);
1243         }
1244     }
1245     if (MASTER(cr) && !(*bConverged))
1246     {
1247         /* Note that the energies and virial are incorrect when not converged */
1248         if (fplog)
1249         {
1250             fprintf(fplog,
1251                     "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1252                     gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1253         }
1254         fprintf(stderr,
1255                 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1256                 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1257     }
1258
1259     /* Copy back the coordinates and the forces */
1260     memcpy(state->x, pos[Min], nat*sizeof(state->x[0]));
1261     memcpy(f, force[Min], nat*sizeof(f[0]));
1262
1263     return count;
1264 }