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