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