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