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