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