e5b7a3d75fbb022b1a2114419cb3bc320a1ef53f
[alexxy/gromacs.git] / src / gromacs / mdlib / shellfc.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2008, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include <stdlib.h>
42 #include <string.h>
43
44 #include "gromacs/legacyheaders/typedefs.h"
45 #include "gromacs/legacyheaders/types/commrec.h"
46 #include "gromacs/legacyheaders/txtdump.h"
47 #include "gromacs/legacyheaders/force.h"
48 #include "gromacs/legacyheaders/mdrun.h"
49 #include "gromacs/legacyheaders/mdatoms.h"
50 #include "gromacs/legacyheaders/vsite.h"
51 #include "gromacs/legacyheaders/network.h"
52 #include "gromacs/legacyheaders/names.h"
53 #include "gromacs/legacyheaders/constr.h"
54 #include "gromacs/legacyheaders/domdec.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/legacyheaders/shellfc.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/legacyheaders/chargegroup.h"
59 #include "gromacs/legacyheaders/macros.h"
60
61 #include "gromacs/math/vec.h"
62 #include "gromacs/pbcutil/mshift.h"
63 #include "gromacs/pbcutil/pbc.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, dt_2, dt_3, 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                     m1 = atom[n1].m;
191                     m2 = atom[n2].m;
192                 }
193                 tm = dt_1/(m1+m2);
194                 for (m = 0; (m < DIM); m++)
195                 {
196                     x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
197                 }
198                 break;
199             case 3:
200                 n1 = s[i].nucl1;
201                 n2 = s[i].nucl2;
202                 n3 = s[i].nucl3;
203                 if (mass)
204                 {
205                     m1 = mass[n1];
206                     m2 = mass[n2];
207                     m3 = mass[n3];
208                 }
209                 else
210                 {
211                     /* Not the correct masses with FE, but it is just a prediction... */
212                     gmx_mtop_atomnr_to_atom(alook, n1, &atom);
213                     m1 = atom->m;
214                     gmx_mtop_atomnr_to_atom(alook, n2, &atom);
215                     m2 = atom->m;
216                     gmx_mtop_atomnr_to_atom(alook, n3, &atom);
217                     m3 = atom->m;
218                 }
219                 tm = dt_1/(m1+m2+m3);
220                 for (m = 0; (m < DIM); m++)
221                 {
222                     x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
223                 }
224                 break;
225             default:
226                 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
227         }
228     }
229
230     if (mass == NULL)
231     {
232         gmx_mtop_atomlookup_destroy(alook);
233     }
234 }
235
236 gmx_shellfc_t init_shell_flexcon(FILE *fplog,
237                                  gmx_mtop_t *mtop, int nflexcon,
238                                  rvec *x)
239 {
240     struct gmx_shellfc       *shfc;
241     t_shell                  *shell;
242     int                      *shell_index = NULL, *at2cg;
243     t_atom                   *atom;
244     int                       n[eptNR], ns, nshell, nsi;
245     int                       i, j, nmol, type, mb, mt, a_offset, cg, mol, ftype, nra;
246     real                      qS, alpha;
247     int                       aS, aN = 0; /* Shell and nucleus */
248     int                       bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
249 #define NBT asize(bondtypes)
250     t_iatom                  *ia;
251     gmx_mtop_atomloop_block_t aloopb;
252     gmx_mtop_atomloop_all_t   aloop;
253     gmx_ffparams_t           *ffparams;
254     gmx_molblock_t           *molb;
255     gmx_moltype_t            *molt;
256     t_block                  *cgs;
257
258     /* Count number of shells, and find their indices */
259     for (i = 0; (i < eptNR); i++)
260     {
261         n[i] = 0;
262     }
263
264     aloopb = gmx_mtop_atomloop_block_init(mtop);
265     while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
266     {
267         n[atom->ptype] += nmol;
268     }
269
270     if (fplog)
271     {
272         /* Print the number of each particle type */
273         for (i = 0; (i < eptNR); i++)
274         {
275             if (n[i] != 0)
276             {
277                 fprintf(fplog, "There are: %d %ss\n", n[i], ptype_str[i]);
278             }
279         }
280     }
281
282     nshell = n[eptShell];
283
284     if (nshell == 0 && nflexcon == 0)
285     {
286         /* We're not doing shells or flexible constraints */
287         return NULL;
288     }
289
290     snew(shfc, 1);
291     shfc->nflexcon = nflexcon;
292
293     if (nshell == 0)
294     {
295         return shfc;
296     }
297
298     /* We have shells: fill the shell data structure */
299
300     /* Global system sized array, this should be avoided */
301     snew(shell_index, mtop->natoms);
302
303     aloop  = gmx_mtop_atomloop_all_init(mtop);
304     nshell = 0;
305     while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
306     {
307         if (atom->ptype == eptShell)
308         {
309             shell_index[i] = nshell++;
310         }
311     }
312
313     snew(shell, nshell);
314
315     /* Initiate the shell structures */
316     for (i = 0; (i < nshell); i++)
317     {
318         shell[i].shell = NO_ATID;
319         shell[i].nnucl = 0;
320         shell[i].nucl1 = NO_ATID;
321         shell[i].nucl2 = NO_ATID;
322         shell[i].nucl3 = NO_ATID;
323         /* shell[i].bInterCG=FALSE; */
324         shell[i].k_1   = 0;
325         shell[i].k     = 0;
326     }
327
328     ffparams = &mtop->ffparams;
329
330     /* Now fill the structures */
331     shfc->bInterCG = FALSE;
332     ns             = 0;
333     a_offset       = 0;
334     for (mb = 0; mb < mtop->nmolblock; mb++)
335     {
336         molb = &mtop->molblock[mb];
337         molt = &mtop->moltype[molb->type];
338
339         cgs = &molt->cgs;
340         snew(at2cg, molt->atoms.nr);
341         for (cg = 0; cg < cgs->nr; cg++)
342         {
343             for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
344             {
345                 at2cg[i] = cg;
346             }
347         }
348
349         atom = molt->atoms.atom;
350         for (mol = 0; mol < molb->nmol; mol++)
351         {
352             for (j = 0; (j < NBT); j++)
353             {
354                 ia = molt->ilist[bondtypes[j]].iatoms;
355                 for (i = 0; (i < molt->ilist[bondtypes[j]].nr); )
356                 {
357                     type  = ia[0];
358                     ftype = ffparams->functype[type];
359                     nra   = interaction_function[ftype].nratoms;
360
361                     /* Check whether we have a bond with a shell */
362                     aS = NO_ATID;
363
364                     switch (bondtypes[j])
365                     {
366                         case F_BONDS:
367                         case F_HARMONIC:
368                         case F_CUBICBONDS:
369                         case F_POLARIZATION:
370                         case F_ANHARM_POL:
371                             if (atom[ia[1]].ptype == eptShell)
372                             {
373                                 aS = ia[1];
374                                 aN = ia[2];
375                             }
376                             else if (atom[ia[2]].ptype == eptShell)
377                             {
378                                 aS = ia[2];
379                                 aN = ia[1];
380                             }
381                             break;
382                         case F_WATER_POL:
383                             aN    = ia[4]; /* Dummy */
384                             aS    = ia[5]; /* Shell */
385                             break;
386                         default:
387                             gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
388                     }
389
390                     if (aS != NO_ATID)
391                     {
392                         qS = atom[aS].q;
393
394                         /* Check whether one of the particles is a shell... */
395                         nsi = shell_index[a_offset+aS];
396                         if ((nsi < 0) || (nsi >= nshell))
397                         {
398                             gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
399                                       nsi, nshell, aS);
400                         }
401                         if (shell[nsi].shell == NO_ATID)
402                         {
403                             shell[nsi].shell = a_offset + aS;
404                             ns++;
405                         }
406                         else if (shell[nsi].shell != a_offset+aS)
407                         {
408                             gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
409                         }
410
411                         if      (shell[nsi].nucl1 == NO_ATID)
412                         {
413                             shell[nsi].nucl1 = a_offset + aN;
414                         }
415                         else if (shell[nsi].nucl2 == NO_ATID)
416                         {
417                             shell[nsi].nucl2 = a_offset + aN;
418                         }
419                         else if (shell[nsi].nucl3 == NO_ATID)
420                         {
421                             shell[nsi].nucl3 = a_offset + aN;
422                         }
423                         else
424                         {
425                             if (fplog)
426                             {
427                                 pr_shell(fplog, ns, shell);
428                             }
429                             gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
430                         }
431                         if (at2cg[aS] != at2cg[aN])
432                         {
433                             /* shell[nsi].bInterCG = TRUE; */
434                             shfc->bInterCG = TRUE;
435                         }
436
437                         switch (bondtypes[j])
438                         {
439                             case F_BONDS:
440                             case F_HARMONIC:
441                                 shell[nsi].k    += ffparams->iparams[type].harmonic.krA;
442                                 break;
443                             case F_CUBICBONDS:
444                                 shell[nsi].k    += ffparams->iparams[type].cubic.kb;
445                                 break;
446                             case F_POLARIZATION:
447                             case F_ANHARM_POL:
448                                 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
449                                 {
450                                     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);
451                                 }
452                                 shell[nsi].k    += sqr(qS)*ONE_4PI_EPS0/
453                                     ffparams->iparams[type].polarize.alpha;
454                                 break;
455                             case F_WATER_POL:
456                                 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
457                                 {
458                                     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);
459                                 }
460                                 alpha          = (ffparams->iparams[type].wpol.al_x+
461                                                   ffparams->iparams[type].wpol.al_y+
462                                                   ffparams->iparams[type].wpol.al_z)/3.0;
463                                 shell[nsi].k  += sqr(qS)*ONE_4PI_EPS0/alpha;
464                                 break;
465                             default:
466                                 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
467                         }
468                         shell[nsi].nnucl++;
469                     }
470                     ia += nra+1;
471                     i  += nra+1;
472                 }
473             }
474             a_offset += molt->atoms.nr;
475         }
476         /* Done with this molecule type */
477         sfree(at2cg);
478     }
479
480     /* Verify whether it's all correct */
481     if (ns != nshell)
482     {
483         gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
484     }
485
486     for (i = 0; (i < ns); i++)
487     {
488         shell[i].k_1 = 1.0/shell[i].k;
489     }
490
491     if (debug)
492     {
493         pr_shell(debug, ns, shell);
494     }
495
496
497     shfc->nshell_gl      = ns;
498     shfc->shell_gl       = shell;
499     shfc->shell_index_gl = shell_index;
500
501     shfc->bPredict     = (getenv("GMX_NOPREDICT") == NULL);
502     shfc->bRequireInit = FALSE;
503     if (!shfc->bPredict)
504     {
505         if (fplog)
506         {
507             fprintf(fplog, "\nWill never predict shell positions\n");
508         }
509     }
510     else
511     {
512         shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != NULL);
513         if (shfc->bRequireInit && fplog)
514         {
515             fprintf(fplog, "\nWill always initiate shell positions\n");
516         }
517     }
518
519     if (shfc->bPredict)
520     {
521         if (x)
522         {
523             predict_shells(fplog, x, NULL, 0, shfc->nshell_gl, shfc->shell_gl,
524                            NULL, mtop, TRUE);
525         }
526
527         if (shfc->bInterCG)
528         {
529             if (fplog)
530             {
531                 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");
532             }
533             /* Prediction improves performance, so we should implement either:
534              * 1. communication for the atoms needed for prediction
535              * 2. prediction using the velocities of shells; currently the
536              *    shell velocities are zeroed, it's a bit tricky to keep
537              *    track of the shell displacements and thus the velocity.
538              */
539             shfc->bPredict = FALSE;
540         }
541     }
542
543     return shfc;
544 }
545
546 void make_local_shells(t_commrec *cr, t_mdatoms *md,
547                        struct gmx_shellfc *shfc)
548 {
549     t_shell      *shell;
550     int           a0, a1, *ind, nshell, i;
551     gmx_domdec_t *dd = NULL;
552
553     if (DOMAINDECOMP(cr))
554     {
555         dd = cr->dd;
556         a0 = 0;
557         a1 = dd->nat_home;
558     }
559     else
560     {
561         /* Single node: we need all shells, just copy the pointer */
562         shfc->nshell = shfc->nshell_gl;
563         shfc->shell  = shfc->shell_gl;
564
565         return;
566     }
567
568     ind = shfc->shell_index_gl;
569
570     nshell = 0;
571     shell  = shfc->shell;
572     for (i = a0; i < a1; i++)
573     {
574         if (md->ptype[i] == eptShell)
575         {
576             if (nshell+1 > shfc->shell_nalloc)
577             {
578                 shfc->shell_nalloc = over_alloc_dd(nshell+1);
579                 srenew(shell, shfc->shell_nalloc);
580             }
581             if (dd)
582             {
583                 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
584             }
585             else
586             {
587                 shell[nshell] = shfc->shell_gl[ind[i]];
588             }
589
590             /* With inter-cg shells we can no do shell prediction,
591              * so we do not need the nuclei numbers.
592              */
593             if (!shfc->bInterCG)
594             {
595                 shell[nshell].nucl1   = i + shell[nshell].nucl1 - shell[nshell].shell;
596                 if (shell[nshell].nnucl > 1)
597                 {
598                     shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
599                 }
600                 if (shell[nshell].nnucl > 2)
601                 {
602                     shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
603                 }
604             }
605             shell[nshell].shell = i;
606             nshell++;
607         }
608     }
609
610     shfc->nshell = nshell;
611     shfc->shell  = shell;
612 }
613
614 static void do_1pos(rvec xnew, rvec xold, rvec f, real step)
615 {
616     real xo, yo, zo;
617     real dx, dy, dz;
618
619     xo = xold[XX];
620     yo = xold[YY];
621     zo = xold[ZZ];
622
623     dx = f[XX]*step;
624     dy = f[YY]*step;
625     dz = f[ZZ]*step;
626
627     xnew[XX] = xo+dx;
628     xnew[YY] = yo+dy;
629     xnew[ZZ] = zo+dz;
630 }
631
632 static void do_1pos3(rvec xnew, rvec xold, rvec f, rvec step)
633 {
634     real xo, yo, zo;
635     real dx, dy, dz;
636
637     xo = xold[XX];
638     yo = xold[YY];
639     zo = xold[ZZ];
640
641     dx = f[XX]*step[XX];
642     dy = f[YY]*step[YY];
643     dz = f[ZZ]*step[ZZ];
644
645     xnew[XX] = xo+dx;
646     xnew[YY] = yo+dy;
647     xnew[ZZ] = zo+dz;
648 }
649
650 static void directional_sd(rvec xold[], rvec xnew[], rvec acc_dir[],
651                            int start, int homenr, real step)
652 {
653     int  i;
654
655     for (i = start; i < homenr; i++)
656     {
657         do_1pos(xnew[i], xold[i], acc_dir[i], step);
658     }
659 }
660
661 static void shell_pos_sd(rvec xcur[], rvec xnew[], rvec f[],
662                          int ns, t_shell s[], int count)
663 {
664     const real step_scale_min       = 0.8,
665                step_scale_increment = 0.2,
666                step_scale_max       = 1.2,
667                step_scale_multiple  = (step_scale_max - step_scale_min) / step_scale_increment;
668     int  i, shell, d;
669     real dx, df, k_est;
670 #ifdef PRINT_STEP
671     real step_min, step_max;
672
673     step_min = 1e30;
674     step_max = 0;
675 #endif
676     for (i = 0; (i < ns); i++)
677     {
678         shell = s[i].shell;
679         if (count == 1)
680         {
681             for (d = 0; d < DIM; d++)
682             {
683                 s[i].step[d] = s[i].k_1;
684 #ifdef PRINT_STEP
685                 step_min = min(step_min, s[i].step[d]);
686                 step_max = max(step_max, s[i].step[d]);
687 #endif
688             }
689         }
690         else
691         {
692             for (d = 0; d < DIM; d++)
693             {
694                 dx = xcur[shell][d] - s[i].xold[d];
695                 df =    f[shell][d] - s[i].fold[d];
696                 /* -dx/df gets used to generate an interpolated value, but would
697                  * cause a NaN if df were binary-equal to zero. Values close to
698                  * zero won't cause problems (because of the min() and max()), so
699                  * just testing for binary inequality is OK. */
700                 if (0.0 != df)
701                 {
702                     k_est = -dx/df;
703                     /* Scale the step size by a factor interpolated from
704                      * step_scale_min to step_scale_max, as k_est goes from 0 to
705                      * step_scale_multiple * s[i].step[d] */
706                     s[i].step[d] =
707                         step_scale_min * s[i].step[d] +
708                         step_scale_increment * min(step_scale_multiple * s[i].step[d], max(k_est, 0));
709                 }
710                 else
711                 {
712                     /* Here 0 == df */
713                     if (gmx_numzero(dx)) /* 0 == dx */
714                     {
715                         /* Likely this will never happen, but if it does just
716                          * don't scale the step. */
717                     }
718                     else /* 0 != dx */
719                     {
720                         s[i].step[d] *= step_scale_max;
721                     }
722                 }
723 #ifdef PRINT_STEP
724                 step_min = min(step_min, s[i].step[d]);
725                 step_max = max(step_max, s[i].step[d]);
726 #endif
727             }
728         }
729         copy_rvec(xcur[shell], s[i].xold);
730         copy_rvec(f[shell],   s[i].fold);
731
732         do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
733
734         if (gmx_debug_at)
735         {
736             fprintf(debug, "shell[%d] = %d\n", i, shell);
737             pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
738             pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
739             pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
740             pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
741         }
742     }
743 #ifdef PRINT_STEP
744     printf("step %.3e %.3e\n", step_min, step_max);
745 #endif
746 }
747
748 static void decrease_step_size(int nshell, t_shell s[])
749 {
750     int i;
751
752     for (i = 0; i < nshell; i++)
753     {
754         svmul(0.8, s[i].step, s[i].step);
755     }
756 }
757
758 static void print_epot(FILE *fp, gmx_int64_t mdstep, int count, real epot, real df,
759                        int ndir, real sf_dir)
760 {
761     char buf[22];
762
763     fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
764             gmx_step_str(mdstep, buf), count, epot, df);
765     if (ndir)
766     {
767         fprintf(fp, ", dir. rmsF: %6.2e\n", sqrt(sf_dir/ndir));
768     }
769     else
770     {
771         fprintf(fp, "\n");
772     }
773 }
774
775
776 static real rms_force(t_commrec *cr, rvec f[], int ns, t_shell s[],
777                       int ndir, real *sf_dir, real *Epot)
778 {
779     int    i, shell, ntot;
780     double buf[4];
781
782     buf[0] = *sf_dir;
783     for (i = 0; i < ns; i++)
784     {
785         shell    = s[i].shell;
786         buf[0]  += norm2(f[shell]);
787     }
788     ntot = ns;
789
790     if (PAR(cr))
791     {
792         buf[1] = ntot;
793         buf[2] = *sf_dir;
794         buf[3] = *Epot;
795         gmx_sumd(4, buf, cr);
796         ntot    = (int)(buf[1] + 0.5);
797         *sf_dir = buf[2];
798         *Epot   = buf[3];
799     }
800     ntot += ndir;
801
802     return (ntot ? sqrt(buf[0]/ntot) : 0);
803 }
804
805 static void check_pbc(FILE *fp, rvec x[], int shell)
806 {
807     int m, now;
808
809     now = shell-4;
810     for (m = 0; (m < DIM); m++)
811     {
812         if (fabs(x[shell][m]-x[now][m]) > 0.3)
813         {
814             pr_rvecs(fp, 0, "SHELL-X", x+now, 5);
815             break;
816         }
817     }
818 }
819
820 static void dump_shells(FILE *fp, rvec x[], rvec f[], real ftol, int ns, t_shell s[])
821 {
822     int  i, shell;
823     real ft2, ff2;
824
825     ft2 = sqr(ftol);
826
827     for (i = 0; (i < ns); i++)
828     {
829         shell = s[i].shell;
830         ff2   = iprod(f[shell], f[shell]);
831         if (ff2 > ft2)
832         {
833             fprintf(fp, "SHELL %5d, force %10.5f  %10.5f  %10.5f, |f| %10.5f\n",
834                     shell, f[shell][XX], f[shell][YY], f[shell][ZZ], sqrt(ff2));
835         }
836         check_pbc(fp, x, shell);
837     }
838 }
839
840 static void init_adir(FILE *log, gmx_shellfc_t shfc,
841                       gmx_constr_t constr, t_idef *idef, t_inputrec *ir,
842                       t_commrec *cr, int dd_ac1,
843                       gmx_int64_t step, t_mdatoms *md, int start, int end,
844                       rvec *x_old, rvec *x_init, rvec *x,
845                       rvec *f, rvec *acc_dir,
846                       gmx_bool bMolPBC, matrix box,
847                       real *lambda, real *dvdlambda, t_nrnb *nrnb)
848 {
849     rvec           *xnold, *xnew;
850     double          w_dt;
851     int             gf, ga, gt;
852     real            dt, scale;
853     int             n, d;
854     unsigned short *ptype;
855     rvec            p, dx;
856
857     if (DOMAINDECOMP(cr))
858     {
859         n = dd_ac1;
860     }
861     else
862     {
863         n = end - start;
864     }
865     if (n > shfc->adir_nalloc)
866     {
867         shfc->adir_nalloc = over_alloc_dd(n);
868         srenew(shfc->adir_xnold, shfc->adir_nalloc);
869         srenew(shfc->adir_xnew, shfc->adir_nalloc);
870     }
871     xnold = shfc->adir_xnold;
872     xnew  = shfc->adir_xnew;
873
874     ptype = md->ptype;
875
876     dt = ir->delta_t;
877
878     /* Does NOT work with freeze or acceleration groups (yet) */
879     for (n = start; n < end; n++)
880     {
881         w_dt = md->invmass[n]*dt;
882
883         for (d = 0; d < DIM; d++)
884         {
885             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
886             {
887                 xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
888                 xnew[n-start][d]  = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
889             }
890             else
891             {
892                 xnold[n-start][d] = x[n][d];
893                 xnew[n-start][d]  = x[n][d];
894             }
895         }
896     }
897     constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, 1.0, md,
898               x, xnold-start, NULL, bMolPBC, box,
899               lambda[efptBONDED], &(dvdlambda[efptBONDED]),
900               NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
901     constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, 1.0, md,
902               x, xnew-start, NULL, bMolPBC, box,
903               lambda[efptBONDED], &(dvdlambda[efptBONDED]),
904               NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
905
906     for (n = start; n < end; n++)
907     {
908         for (d = 0; d < DIM; d++)
909         {
910             xnew[n-start][d] =
911                 -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
912                 - f[n][d]*md->invmass[n];
913         }
914         clear_rvec(acc_dir[n]);
915     }
916
917     /* Project the acceleration on the old bond directions */
918     constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, 1.0, md,
919               x_old, xnew-start, acc_dir, bMolPBC, box,
920               lambda[efptBONDED], &(dvdlambda[efptBONDED]),
921               NULL, NULL, nrnb, econqDeriv_FlexCon, FALSE, 0, 0);
922 }
923
924 int relax_shell_flexcon(FILE *fplog, t_commrec *cr, gmx_bool bVerbose,
925                         gmx_int64_t mdstep, t_inputrec *inputrec,
926                         gmx_bool bDoNS, int force_flags,
927                         gmx_localtop_t *top,
928                         gmx_constr_t constr,
929                         gmx_enerdata_t *enerd, t_fcdata *fcd,
930                         t_state *state, rvec f[],
931                         tensor force_vir,
932                         t_mdatoms *md,
933                         t_nrnb *nrnb, gmx_wallcycle_t wcycle,
934                         t_graph *graph,
935                         gmx_groups_t *groups,
936                         struct gmx_shellfc *shfc,
937                         t_forcerec *fr,
938                         gmx_bool bBornRadii,
939                         double t, rvec mu_tot,
940                         gmx_bool *bConverged,
941                         gmx_vsite_t *vsite,
942                         FILE *fp_field)
943 {
944     int        nshell;
945     t_shell   *shell;
946     t_idef    *idef;
947     rvec      *pos[2], *force[2], *acc_dir = NULL, *x_old = NULL;
948     real       Epot[2], df[2];
949     rvec       dx;
950     real       sf_dir, invdt;
951     real       ftol, xiH, xiS, 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, g, 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 = 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 }