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