Merge branch release-5-0 into release-5-1
[alexxy/gromacs.git] / src / gromacs / mdlib / shellfc.cpp
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 #include "gmxpre.h"
38
39 #include "gromacs/legacyheaders/shellfc.h"
40
41 #include <stdlib.h>
42 #include <string.h>
43
44 #include <algorithm>
45
46 #include "gromacs/domdec/domdec.h"
47 #include "gromacs/legacyheaders/chargegroup.h"
48 #include "gromacs/legacyheaders/constr.h"
49 #include "gromacs/legacyheaders/force.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/mdatoms.h"
52 #include "gromacs/legacyheaders/mdrun.h"
53 #include "gromacs/legacyheaders/names.h"
54 #include "gromacs/legacyheaders/network.h"
55 #include "gromacs/legacyheaders/txtdump.h"
56 #include "gromacs/legacyheaders/typedefs.h"
57 #include "gromacs/legacyheaders/vsite.h"
58 #include "gromacs/legacyheaders/types/commrec.h"
59 #include "gromacs/math/units.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/pbcutil/mshift.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/smalloc.h"
66
67 typedef struct {
68     int     nnucl;
69     atom_id shell;               /* The shell id                                */
70     atom_id nucl1, nucl2, nucl3; /* The nuclei connected to the shell   */
71     /* gmx_bool    bInterCG; */       /* Coupled to nuclei outside cg?        */
72     real    k;                   /* force constant                      */
73     real    k_1;                 /* 1 over force constant               */
74     rvec    xold;
75     rvec    fold;
76     rvec    step;
77 } t_shell;
78
79 typedef struct gmx_shellfc {
80     int         nshell_gl;      /* The number of shells in the system       */
81     t_shell    *shell_gl;       /* All the shells (for DD only)             */
82     int        *shell_index_gl; /* Global shell index (for DD only)         */
83     gmx_bool    bInterCG;       /* Are there inter charge-group shells?     */
84     int         nshell;         /* The number of local shells               */
85     t_shell    *shell;          /* The local shells                         */
86     int         shell_nalloc;   /* The allocation size of shell             */
87     gmx_bool    bPredict;       /* Predict shell positions                  */
88     gmx_bool    bRequireInit;   /* Require initialization of shell positions  */
89     int         nflexcon;       /* The number of flexible constraints       */
90     rvec       *x[2];           /* Array for iterative minimization         */
91     rvec       *f[2];           /* Array for iterative minimization         */
92     int         x_nalloc;       /* The allocation size of x and f           */
93     rvec       *acc_dir;        /* Acceleration direction for flexcon       */
94     rvec       *x_old;          /* Old coordinates for flexcon              */
95     int         flex_nalloc;    /* The allocation size of acc_dir and x_old */
96     rvec       *adir_xnold;     /* Work space for init_adir                 */
97     rvec       *adir_xnew;      /* Work space for init_adir                 */
98     int         adir_nalloc;    /* Work space for init_adir                 */
99 } t_gmx_shellfc;
100
101
102 static void pr_shell(FILE *fplog, int ns, t_shell s[])
103 {
104     int i;
105
106     fprintf(fplog, "SHELL DATA\n");
107     fprintf(fplog, "%5s  %8s  %5s  %5s  %5s\n",
108             "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
109     for (i = 0; (i < ns); i++)
110     {
111         fprintf(fplog, "%5d  %8.3f  %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
112         if (s[i].nnucl == 2)
113         {
114             fprintf(fplog, "  %5d\n", s[i].nucl2);
115         }
116         else if (s[i].nnucl == 3)
117         {
118             fprintf(fplog, "  %5d  %5d\n", s[i].nucl2, s[i].nucl3);
119         }
120         else
121         {
122             fprintf(fplog, "\n");
123         }
124     }
125 }
126
127 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
128                            int ns, t_shell s[],
129                            real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
130 {
131     int                   i, m, s1, n1, n2, n3;
132     real                  dt_1, fudge, tm, m1, m2, m3;
133     rvec                 *ptr;
134     gmx_mtop_atomlookup_t alook = NULL;
135     t_atom               *atom;
136
137     if (mass == NULL)
138     {
139         alook = gmx_mtop_atomlookup_init(mtop);
140     }
141
142     /* We introduce a fudge factor for performance reasons: with this choice
143      * the initial force on the shells is about a factor of two lower than
144      * without
145      */
146     fudge = 1.0;
147
148     if (bInit)
149     {
150         if (fplog)
151         {
152             fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
153         }
154         ptr  = x;
155         dt_1 = 1;
156     }
157     else
158     {
159         ptr  = v;
160         dt_1 = fudge*dt;
161     }
162
163     for (i = 0; (i < ns); i++)
164     {
165         s1 = s[i].shell;
166         if (bInit)
167         {
168             clear_rvec(x[s1]);
169         }
170         switch (s[i].nnucl)
171         {
172             case 1:
173                 n1 = s[i].nucl1;
174                 for (m = 0; (m < DIM); m++)
175                 {
176                     x[s1][m] += ptr[n1][m]*dt_1;
177                 }
178                 break;
179             case 2:
180                 n1 = s[i].nucl1;
181                 n2 = s[i].nucl2;
182                 if (mass)
183                 {
184                     m1 = mass[n1];
185                     m2 = mass[n2];
186                 }
187                 else
188                 {
189                     /* Not the correct masses with FE, but it is just a prediction... */
190                     gmx_mtop_atomnr_to_atom(alook, n1, &atom);
191                     m1 = atom->m;
192                     gmx_mtop_atomnr_to_atom(alook, n2, &atom);
193                     m2 = atom->m;
194                 }
195                 tm = dt_1/(m1+m2);
196                 for (m = 0; (m < DIM); m++)
197                 {
198                     x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
199                 }
200                 break;
201             case 3:
202                 n1 = s[i].nucl1;
203                 n2 = s[i].nucl2;
204                 n3 = s[i].nucl3;
205                 if (mass)
206                 {
207                     m1 = mass[n1];
208                     m2 = mass[n2];
209                     m3 = mass[n3];
210                 }
211                 else
212                 {
213                     /* Not the correct masses with FE, but it is just a prediction... */
214                     gmx_mtop_atomnr_to_atom(alook, n1, &atom);
215                     m1 = atom->m;
216                     gmx_mtop_atomnr_to_atom(alook, n2, &atom);
217                     m2 = atom->m;
218                     gmx_mtop_atomnr_to_atom(alook, n3, &atom);
219                     m3 = atom->m;
220                 }
221                 tm = dt_1/(m1+m2+m3);
222                 for (m = 0; (m < DIM); m++)
223                 {
224                     x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
225                 }
226                 break;
227             default:
228                 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
229         }
230     }
231
232     if (mass == NULL)
233     {
234         gmx_mtop_atomlookup_destroy(alook);
235     }
236 }
237
238 gmx_shellfc_t init_shell_flexcon(FILE *fplog,
239                                  gmx_mtop_t *mtop, int nflexcon,
240                                  rvec *x)
241 {
242     struct gmx_shellfc       *shfc;
243     t_shell                  *shell;
244     int                      *shell_index = NULL, *at2cg;
245     t_atom                   *atom;
246     int                       n[eptNR], ns, nshell, nsi;
247     int                       i, j, nmol, type, mb, a_offset, cg, mol, ftype, nra;
248     real                      qS, alpha;
249     int                       aS, aN = 0; /* Shell and nucleus */
250     int                       bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
251 #define NBT asize(bondtypes)
252     t_iatom                  *ia;
253     gmx_mtop_atomloop_block_t aloopb;
254     gmx_mtop_atomloop_all_t   aloop;
255     gmx_ffparams_t           *ffparams;
256     gmx_molblock_t           *molb;
257     gmx_moltype_t            *molt;
258     t_block                  *cgs;
259
260     /* Count number of shells, and find their indices */
261     for (i = 0; (i < eptNR); i++)
262     {
263         n[i] = 0;
264     }
265
266     aloopb = gmx_mtop_atomloop_block_init(mtop);
267     while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
268     {
269         n[atom->ptype] += nmol;
270     }
271
272     if (fplog)
273     {
274         /* Print the number of each particle type */
275         for (i = 0; (i < eptNR); i++)
276         {
277             if (n[i] != 0)
278             {
279                 fprintf(fplog, "There are: %d %ss\n", n[i], ptype_str[i]);
280             }
281         }
282     }
283
284     nshell = n[eptShell];
285
286     if (nshell == 0 && nflexcon == 0)
287     {
288         /* We're not doing shells or flexible constraints */
289         return NULL;
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             /* Prediction improves performance, so we should implement either:
536              * 1. communication for the atoms needed for prediction
537              * 2. prediction using the velocities of shells; currently the
538              *    shell velocities are zeroed, it's a bit tricky to keep
539              *    track of the shell displacements and thus the velocity.
540              */
541             shfc->bPredict = FALSE;
542         }
543     }
544
545     return shfc;
546 }
547
548 void make_local_shells(t_commrec *cr, t_mdatoms *md,
549                        struct gmx_shellfc *shfc)
550 {
551     t_shell      *shell;
552     int           a0, a1, *ind, nshell, i;
553     gmx_domdec_t *dd = NULL;
554
555     if (DOMAINDECOMP(cr))
556     {
557         dd = cr->dd;
558         a0 = 0;
559         a1 = dd->nat_home;
560     }
561     else
562     {
563         /* Single node: we need all shells, just copy the pointer */
564         shfc->nshell = shfc->nshell_gl;
565         shfc->shell  = shfc->shell_gl;
566
567         return;
568     }
569
570     ind = shfc->shell_index_gl;
571
572     nshell = 0;
573     shell  = shfc->shell;
574     for (i = a0; i < a1; i++)
575     {
576         if (md->ptype[i] == eptShell)
577         {
578             if (nshell+1 > shfc->shell_nalloc)
579             {
580                 shfc->shell_nalloc = over_alloc_dd(nshell+1);
581                 srenew(shell, shfc->shell_nalloc);
582             }
583             if (dd)
584             {
585                 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
586             }
587             else
588             {
589                 shell[nshell] = shfc->shell_gl[ind[i]];
590             }
591
592             /* With inter-cg shells we can no do shell prediction,
593              * so we do not need the nuclei numbers.
594              */
595             if (!shfc->bInterCG)
596             {
597                 shell[nshell].nucl1   = i + shell[nshell].nucl1 - shell[nshell].shell;
598                 if (shell[nshell].nnucl > 1)
599                 {
600                     shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
601                 }
602                 if (shell[nshell].nnucl > 2)
603                 {
604                     shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
605                 }
606             }
607             shell[nshell].shell = i;
608             nshell++;
609         }
610     }
611
612     shfc->nshell = nshell;
613     shfc->shell  = shell;
614 }
615
616 static void do_1pos(rvec xnew, rvec xold, rvec f, real step)
617 {
618     real xo, yo, zo;
619     real dx, dy, dz;
620
621     xo = xold[XX];
622     yo = xold[YY];
623     zo = xold[ZZ];
624
625     dx = f[XX]*step;
626     dy = f[YY]*step;
627     dz = f[ZZ]*step;
628
629     xnew[XX] = xo+dx;
630     xnew[YY] = yo+dy;
631     xnew[ZZ] = zo+dz;
632 }
633
634 static void do_1pos3(rvec xnew, rvec xold, rvec f, rvec step)
635 {
636     real xo, yo, zo;
637     real dx, dy, dz;
638
639     xo = xold[XX];
640     yo = xold[YY];
641     zo = xold[ZZ];
642
643     dx = f[XX]*step[XX];
644     dy = f[YY]*step[YY];
645     dz = f[ZZ]*step[ZZ];
646
647     xnew[XX] = xo+dx;
648     xnew[YY] = yo+dy;
649     xnew[ZZ] = zo+dz;
650 }
651
652 static void directional_sd(rvec xold[], rvec xnew[], rvec acc_dir[],
653                            int start, int homenr, real step)
654 {
655     int  i;
656
657     for (i = start; i < homenr; i++)
658     {
659         do_1pos(xnew[i], xold[i], acc_dir[i], step);
660     }
661 }
662
663 static void shell_pos_sd(rvec xcur[], rvec xnew[], rvec f[],
664                          int ns, t_shell s[], int count)
665 {
666     const real step_scale_min       = 0.8,
667                step_scale_increment = 0.2,
668                step_scale_max       = 1.2,
669                step_scale_multiple  = (step_scale_max - step_scale_min) / step_scale_increment;
670     int        i, shell, d;
671     real       dx, df, k_est;
672     const real zero = 0;
673 #ifdef PRINT_STEP
674     real       step_min, step_max;
675
676     step_min = 1e30;
677     step_max = 0;
678 #endif
679     for (i = 0; (i < ns); i++)
680     {
681         shell = s[i].shell;
682         if (count == 1)
683         {
684             for (d = 0; d < DIM; d++)
685             {
686                 s[i].step[d] = s[i].k_1;
687 #ifdef PRINT_STEP
688                 step_min = std::min(step_min, s[i].step[d]);
689                 step_max = std::max(step_max, s[i].step[d]);
690 #endif
691             }
692         }
693         else
694         {
695             for (d = 0; d < DIM; d++)
696             {
697                 dx = xcur[shell][d] - s[i].xold[d];
698                 df =    f[shell][d] - s[i].fold[d];
699                 /* -dx/df gets used to generate an interpolated value, but would
700                  * cause a NaN if df were binary-equal to zero. Values close to
701                  * zero won't cause problems (because of the min() and max()), so
702                  * just testing for binary inequality is OK. */
703                 if (0.0 != df)
704                 {
705                     k_est = -dx/df;
706                     /* Scale the step size by a factor interpolated from
707                      * step_scale_min to step_scale_max, as k_est goes from 0 to
708                      * step_scale_multiple * s[i].step[d] */
709                     s[i].step[d] =
710                         step_scale_min * s[i].step[d] +
711                         step_scale_increment * std::min(step_scale_multiple * s[i].step[d], std::max(k_est, zero));
712                 }
713                 else
714                 {
715                     /* Here 0 == df */
716                     if (gmx_numzero(dx)) /* 0 == dx */
717                     {
718                         /* Likely this will never happen, but if it does just
719                          * don't scale the step. */
720                     }
721                     else /* 0 != dx */
722                     {
723                         s[i].step[d] *= step_scale_max;
724                     }
725                 }
726 #ifdef PRINT_STEP
727                 step_min = std::min(step_min, s[i].step[d]);
728                 step_max = std::max(step_max, s[i].step[d]);
729 #endif
730             }
731         }
732         copy_rvec(xcur[shell], s[i].xold);
733         copy_rvec(f[shell],   s[i].fold);
734
735         do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
736
737         if (gmx_debug_at)
738         {
739             fprintf(debug, "shell[%d] = %d\n", i, shell);
740             pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
741             pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
742             pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
743             pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
744         }
745     }
746 #ifdef PRINT_STEP
747     printf("step %.3e %.3e\n", step_min, step_max);
748 #endif
749 }
750
751 static void decrease_step_size(int nshell, t_shell s[])
752 {
753     int i;
754
755     for (i = 0; i < nshell; i++)
756     {
757         svmul(0.8, s[i].step, s[i].step);
758     }
759 }
760
761 static void print_epot(FILE *fp, gmx_int64_t mdstep, int count, real epot, real df,
762                        int ndir, real sf_dir)
763 {
764     char buf[22];
765
766     fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
767             gmx_step_str(mdstep, buf), count, epot, df);
768     if (ndir)
769     {
770         fprintf(fp, ", dir. rmsF: %6.2e\n", sqrt(sf_dir/ndir));
771     }
772     else
773     {
774         fprintf(fp, "\n");
775     }
776 }
777
778
779 static real rms_force(t_commrec *cr, rvec f[], int ns, t_shell s[],
780                       int ndir, real *sf_dir, real *Epot)
781 {
782     int    i, shell, ntot;
783     double buf[4];
784
785     buf[0] = *sf_dir;
786     for (i = 0; i < ns; i++)
787     {
788         shell    = s[i].shell;
789         buf[0]  += norm2(f[shell]);
790     }
791     ntot = ns;
792
793     if (PAR(cr))
794     {
795         buf[1] = ntot;
796         buf[2] = *sf_dir;
797         buf[3] = *Epot;
798         gmx_sumd(4, buf, cr);
799         ntot    = (int)(buf[1] + 0.5);
800         *sf_dir = buf[2];
801         *Epot   = buf[3];
802     }
803     ntot += ndir;
804
805     return (ntot ? sqrt(buf[0]/ntot) : 0);
806 }
807
808 static void check_pbc(FILE *fp, rvec x[], int shell)
809 {
810     int m, now;
811
812     now = shell-4;
813     for (m = 0; (m < DIM); m++)
814     {
815         if (fabs(x[shell][m]-x[now][m]) > 0.3)
816         {
817             pr_rvecs(fp, 0, "SHELL-X", x+now, 5);
818             break;
819         }
820     }
821 }
822
823 static void dump_shells(FILE *fp, rvec x[], rvec f[], real ftol, int ns, t_shell s[])
824 {
825     int  i, shell;
826     real ft2, ff2;
827
828     ft2 = sqr(ftol);
829
830     for (i = 0; (i < ns); i++)
831     {
832         shell = s[i].shell;
833         ff2   = iprod(f[shell], f[shell]);
834         if (ff2 > ft2)
835         {
836             fprintf(fp, "SHELL %5d, force %10.5f  %10.5f  %10.5f, |f| %10.5f\n",
837                     shell, f[shell][XX], f[shell][YY], f[shell][ZZ], sqrt(ff2));
838         }
839         check_pbc(fp, x, shell);
840     }
841 }
842
843 static void init_adir(FILE *log, gmx_shellfc_t shfc,
844                       gmx_constr_t constr, t_idef *idef, t_inputrec *ir,
845                       t_commrec *cr, int dd_ac1,
846                       gmx_int64_t step, t_mdatoms *md, int start, int end,
847                       rvec *x_old, rvec *x_init, rvec *x,
848                       rvec *f, rvec *acc_dir,
849                       gmx_bool bMolPBC, matrix box,
850                       real *lambda, real *dvdlambda, t_nrnb *nrnb)
851 {
852     rvec           *xnold, *xnew;
853     double          w_dt;
854     real            dt;
855     int             n, d;
856     unsigned short *ptype;
857
858     if (DOMAINDECOMP(cr))
859     {
860         n = dd_ac1;
861     }
862     else
863     {
864         n = end - start;
865     }
866     if (n > shfc->adir_nalloc)
867     {
868         shfc->adir_nalloc = over_alloc_dd(n);
869         srenew(shfc->adir_xnold, shfc->adir_nalloc);
870         srenew(shfc->adir_xnew, shfc->adir_nalloc);
871     }
872     xnold = shfc->adir_xnold;
873     xnew  = shfc->adir_xnew;
874
875     ptype = md->ptype;
876
877     dt = ir->delta_t;
878
879     /* Does NOT work with freeze or acceleration groups (yet) */
880     for (n = start; n < end; n++)
881     {
882         w_dt = md->invmass[n]*dt;
883
884         for (d = 0; d < DIM; d++)
885         {
886             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
887             {
888                 xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
889                 xnew[n-start][d]  = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
890             }
891             else
892             {
893                 xnold[n-start][d] = x[n][d];
894                 xnew[n-start][d]  = x[n][d];
895             }
896         }
897     }
898     constrain(log, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
899               x, xnold-start, NULL, bMolPBC, box,
900               lambda[efptBONDED], &(dvdlambda[efptBONDED]),
901               NULL, NULL, nrnb, econqCoord);
902     constrain(log, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
903               x, xnew-start, NULL, bMolPBC, box,
904               lambda[efptBONDED], &(dvdlambda[efptBONDED]),
905               NULL, NULL, nrnb, econqCoord);
906
907     for (n = start; n < end; n++)
908     {
909         for (d = 0; d < DIM; d++)
910         {
911             xnew[n-start][d] =
912                 -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
913                 - f[n][d]*md->invmass[n];
914         }
915         clear_rvec(acc_dir[n]);
916     }
917
918     /* Project the acceleration on the old bond directions */
919     constrain(log, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
920               x_old, xnew-start, acc_dir, bMolPBC, box,
921               lambda[efptBONDED], &(dvdlambda[efptBONDED]),
922               NULL, NULL, nrnb, econqDeriv_FlexCon);
923 }
924
925 int relax_shell_flexcon(FILE *fplog, t_commrec *cr, gmx_bool bVerbose,
926                         gmx_int64_t mdstep, t_inputrec *inputrec,
927                         gmx_bool bDoNS, int force_flags,
928                         gmx_localtop_t *top,
929                         gmx_constr_t constr,
930                         gmx_enerdata_t *enerd, t_fcdata *fcd,
931                         t_state *state, rvec f[],
932                         tensor force_vir,
933                         t_mdatoms *md,
934                         t_nrnb *nrnb, gmx_wallcycle_t wcycle,
935                         t_graph *graph,
936                         gmx_groups_t *groups,
937                         struct gmx_shellfc *shfc,
938                         t_forcerec *fr,
939                         gmx_bool bBornRadii,
940                         double t, rvec mu_tot,
941                         gmx_bool *bConverged,
942                         gmx_vsite_t *vsite,
943                         FILE *fp_field)
944 {
945     int        nshell;
946     t_shell   *shell;
947     t_idef    *idef;
948     rvec      *pos[2], *force[2], *acc_dir = NULL, *x_old = NULL;
949     real       Epot[2], df[2];
950     real       sf_dir, invdt;
951     real       ftol, dum = 0;
952     char       sbuf[22];
953     gmx_bool   bCont, bInit;
954     int        nat, dd_ac0, dd_ac1 = 0, i;
955     int        start = 0, homenr = md->homenr, end = start+homenr, cg0, cg1;
956     int        nflexcon, number_steps, d, Min = 0, count = 0;
957 #define  Try (1-Min)             /* At start Try = 1 */
958
959     bCont        = (mdstep == inputrec->init_step) && inputrec->bContinuation;
960     bInit        = (mdstep == inputrec->init_step) || shfc->bRequireInit;
961     ftol         = inputrec->em_tol;
962     number_steps = inputrec->niter;
963     nshell       = shfc->nshell;
964     shell        = shfc->shell;
965     nflexcon     = shfc->nflexcon;
966
967     idef = &top->idef;
968
969     if (DOMAINDECOMP(cr))
970     {
971         nat = dd_natoms_vsite(cr->dd);
972         if (nflexcon > 0)
973         {
974             dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
975             nat = std::max(nat, dd_ac1);
976         }
977     }
978     else
979     {
980         nat = state->natoms;
981     }
982
983     if (nat > shfc->x_nalloc)
984     {
985         /* Allocate local arrays */
986         shfc->x_nalloc = over_alloc_dd(nat);
987         for (i = 0; (i < 2); i++)
988         {
989             srenew(shfc->x[i], shfc->x_nalloc);
990             srenew(shfc->f[i], shfc->x_nalloc);
991         }
992     }
993     for (i = 0; (i < 2); i++)
994     {
995         pos[i]   = shfc->x[i];
996         force[i] = shfc->f[i];
997     }
998
999     if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1000     {
1001         /* This is the only time where the coordinates are used
1002          * before do_force is called, which normally puts all
1003          * charge groups in the box.
1004          */
1005         if (inputrec->cutoff_scheme == ecutsVERLET)
1006         {
1007             put_atoms_in_box_omp(fr->ePBC, state->box, md->homenr, state->x);
1008         }
1009         else
1010         {
1011             cg0 = 0;
1012             cg1 = top->cgs.nr;
1013             put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1014                                      &(top->cgs), state->x, fr->cg_cm);
1015         }
1016
1017         if (graph)
1018         {
1019             mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
1020         }
1021     }
1022
1023     /* After this all coordinate arrays will contain whole charge groups */
1024     if (graph)
1025     {
1026         shift_self(graph, state->box, state->x);
1027     }
1028
1029     if (nflexcon)
1030     {
1031         if (nat > shfc->flex_nalloc)
1032         {
1033             shfc->flex_nalloc = over_alloc_dd(nat);
1034             srenew(shfc->acc_dir, shfc->flex_nalloc);
1035             srenew(shfc->x_old, shfc->flex_nalloc);
1036         }
1037         acc_dir = shfc->acc_dir;
1038         x_old   = shfc->x_old;
1039         for (i = 0; i < homenr; i++)
1040         {
1041             for (d = 0; d < DIM; d++)
1042             {
1043                 shfc->x_old[i][d] =
1044                     state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
1045             }
1046         }
1047     }
1048
1049     /* Do a prediction of the shell positions */
1050     if (shfc->bPredict && !bCont)
1051     {
1052         predict_shells(fplog, state->x, state->v, inputrec->delta_t, nshell, shell,
1053                        md->massT, NULL, bInit);
1054     }
1055
1056     /* do_force expected the charge groups to be in the box */
1057     if (graph)
1058     {
1059         unshift_self(graph, state->box, state->x);
1060     }
1061
1062     /* Calculate the forces first time around */
1063     if (gmx_debug_at)
1064     {
1065         pr_rvecs(debug, 0, "x b4 do_force", state->x + start, homenr);
1066     }
1067     do_force(fplog, cr, inputrec, mdstep, nrnb, wcycle, top, groups,
1068              state->box, state->x, &state->hist,
1069              force[Min], force_vir, md, enerd, fcd,
1070              state->lambda, graph,
1071              fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1072              (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
1073
1074     sf_dir = 0;
1075     if (nflexcon)
1076     {
1077         init_adir(fplog, shfc,
1078                   constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1079                   shfc->x_old-start, state->x, state->x, force[Min],
1080                   shfc->acc_dir-start,
1081                   fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1082
1083         for (i = start; i < end; i++)
1084         {
1085             sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
1086         }
1087     }
1088
1089     Epot[Min] = enerd->term[F_EPOT];
1090
1091     df[Min] = rms_force(cr, shfc->f[Min], nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1092     df[Try] = 0;
1093     if (debug)
1094     {
1095         fprintf(debug, "df = %g  %g\n", df[Min], df[Try]);
1096     }
1097
1098     if (gmx_debug_at)
1099     {
1100         pr_rvecs(debug, 0, "force0", force[Min], md->nr);
1101     }
1102
1103     if (nshell+nflexcon > 0)
1104     {
1105         /* Copy x to pos[Min] & pos[Try]: during minimization only the
1106          * shell positions are updated, therefore the other particles must
1107          * be set here.
1108          */
1109         memcpy(pos[Min], state->x, nat*sizeof(state->x[0]));
1110         memcpy(pos[Try], state->x, nat*sizeof(state->x[0]));
1111     }
1112
1113     if (bVerbose && MASTER(cr))
1114     {
1115         print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1116     }
1117
1118     if (debug)
1119     {
1120         fprintf(debug, "%17s: %14.10e\n",
1121                 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1122         fprintf(debug, "%17s: %14.10e\n",
1123                 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1124         fprintf(debug, "%17s: %14.10e\n",
1125                 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1126         fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1127     }
1128
1129     /* First check whether we should do shells, or whether the force is
1130      * low enough even without minimization.
1131      */
1132     *bConverged = (df[Min] < ftol);
1133
1134     for (count = 1; (!(*bConverged) && (count < number_steps)); count++)
1135     {
1136         if (vsite)
1137         {
1138             construct_vsites(vsite, pos[Min], inputrec->delta_t, state->v,
1139                              idef->iparams, idef->il,
1140                              fr->ePBC, fr->bMolPBC, cr, state->box);
1141         }
1142
1143         if (nflexcon)
1144         {
1145             init_adir(fplog, shfc,
1146                       constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1147                       x_old-start, state->x, pos[Min], force[Min], acc_dir-start,
1148                       fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1149
1150             directional_sd(pos[Min], pos[Try], acc_dir-start, start, end,
1151                            fr->fc_stepsize);
1152         }
1153
1154         /* New positions, Steepest descent */
1155         shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1156
1157         /* do_force expected the charge groups to be in the box */
1158         if (graph)
1159         {
1160             unshift_self(graph, state->box, pos[Try]);
1161         }
1162
1163         if (gmx_debug_at)
1164         {
1165             pr_rvecs(debug, 0, "RELAX: pos[Min]  ", pos[Min] + start, homenr);
1166             pr_rvecs(debug, 0, "RELAX: pos[Try]  ", pos[Try] + start, homenr);
1167         }
1168         /* Try the new positions */
1169         do_force(fplog, cr, inputrec, 1, nrnb, wcycle,
1170                  top, groups, state->box, pos[Try], &state->hist,
1171                  force[Try], force_vir,
1172                  md, enerd, fcd, state->lambda, graph,
1173                  fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1174                  force_flags);
1175
1176         if (gmx_debug_at)
1177         {
1178             pr_rvecs(debug, 0, "RELAX: force[Min]", force[Min] + start, homenr);
1179             pr_rvecs(debug, 0, "RELAX: force[Try]", force[Try] + start, homenr);
1180         }
1181         sf_dir = 0;
1182         if (nflexcon)
1183         {
1184             init_adir(fplog, shfc,
1185                       constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1186                       x_old-start, state->x, pos[Try], force[Try], acc_dir-start,
1187                       fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1188
1189             for (i = start; i < end; i++)
1190             {
1191                 sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
1192             }
1193         }
1194
1195         Epot[Try] = enerd->term[F_EPOT];
1196
1197         df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1198
1199         if (debug)
1200         {
1201             fprintf(debug, "df = %g  %g\n", df[Min], df[Try]);
1202         }
1203
1204         if (debug)
1205         {
1206             if (gmx_debug_at)
1207             {
1208                 pr_rvecs(debug, 0, "F na do_force", force[Try] + start, homenr);
1209             }
1210             if (gmx_debug_at)
1211             {
1212                 fprintf(debug, "SHELL ITER %d\n", count);
1213                 dump_shells(debug, pos[Try], force[Try], ftol, nshell, shell);
1214             }
1215         }
1216
1217         if (bVerbose && MASTER(cr))
1218         {
1219             print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1220         }
1221
1222         *bConverged = (df[Try] < ftol);
1223
1224         if ((df[Try] < df[Min]))
1225         {
1226             if (debug)
1227             {
1228                 fprintf(debug, "Swapping Min and Try\n");
1229             }
1230             if (nflexcon)
1231             {
1232                 /* Correct the velocities for the flexible constraints */
1233                 invdt = 1/inputrec->delta_t;
1234                 for (i = start; i < end; i++)
1235                 {
1236                     for (d = 0; d < DIM; d++)
1237                     {
1238                         state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1239                     }
1240                 }
1241             }
1242             Min  = Try;
1243         }
1244         else
1245         {
1246             decrease_step_size(nshell, shell);
1247         }
1248     }
1249     if (MASTER(cr) && !(*bConverged))
1250     {
1251         /* Note that the energies and virial are incorrect when not converged */
1252         if (fplog)
1253         {
1254             fprintf(fplog,
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         fprintf(stderr,
1259                 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1260                 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1261     }
1262
1263     /* Copy back the coordinates and the forces */
1264     memcpy(state->x, pos[Min], nat*sizeof(state->x[0]));
1265     memcpy(f, force[Min], nat*sizeof(f[0]));
1266
1267     return count;
1268 }