Propose a post-submit matrix
[gromacs.git] / src / gromacs / mdlib / shellfc.cpp
blobc62ae769c199ccc731b53c417044eb923535625a
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2008, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include "shellfc.h"
41 #include <stdlib.h>
42 #include <string.h>
44 #include <cstdint>
46 #include <algorithm>
47 #include <array>
49 #include "gromacs/domdec/domdec.h"
50 #include "gromacs/domdec/domdec_struct.h"
51 #include "gromacs/gmxlib/chargegroup.h"
52 #include "gromacs/gmxlib/network.h"
53 #include "gromacs/math/functions.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/math/vecdump.h"
57 #include "gromacs/mdlib/constr.h"
58 #include "gromacs/mdlib/force.h"
59 #include "gromacs/mdlib/mdrun.h"
60 #include "gromacs/mdlib/sim_util.h"
61 #include "gromacs/mdlib/vsite.h"
62 #include "gromacs/mdtypes/commrec.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/pbcutil/mshift.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/topology/mtop_util.h"
68 #include "gromacs/utility/arraysize.h"
69 #include "gromacs/utility/cstringutil.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/smalloc.h"
73 typedef struct {
74 int nnucl;
75 int shell; /* The shell id */
76 int nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
77 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
78 real k; /* force constant */
79 real k_1; /* 1 over force constant */
80 rvec xold;
81 rvec fold;
82 rvec step;
83 } t_shell;
85 struct gmx_shellfc_t {
86 int nshell_gl; /* The number of shells in the system */
87 t_shell *shell_gl; /* All the shells (for DD only) */
88 int *shell_index_gl; /* Global shell index (for DD only) */
89 gmx_bool bInterCG; /* Are there inter charge-group shells? */
90 int nshell; /* The number of local shells */
91 t_shell *shell; /* The local shells */
92 int shell_nalloc; /* The allocation size of shell */
93 gmx_bool bPredict; /* Predict shell positions */
94 gmx_bool bRequireInit; /* Require initialization of shell positions */
95 int nflexcon; /* The number of flexible constraints */
96 rvec *x[2]; /* Array for iterative minimization */
97 rvec *f[2]; /* Array for iterative minimization */
98 int x_nalloc; /* The allocation size of x and f */
99 rvec *acc_dir; /* Acceleration direction for flexcon */
100 rvec *x_old; /* Old coordinates for flexcon */
101 int flex_nalloc; /* The allocation size of acc_dir and x_old */
102 rvec *adir_xnold; /* Work space for init_adir */
103 rvec *adir_xnew; /* Work space for init_adir */
104 int adir_nalloc; /* Work space for init_adir */
105 std::int64_t numForceEvaluations; /* Total number of force evaluations */
106 int numConvergedIterations; /* Total number of iterations that converged */
110 static void pr_shell(FILE *fplog, int ns, t_shell s[])
112 int i;
114 fprintf(fplog, "SHELL DATA\n");
115 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
116 "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
117 for (i = 0; (i < ns); i++)
119 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
120 if (s[i].nnucl == 2)
122 fprintf(fplog, " %5d\n", s[i].nucl2);
124 else if (s[i].nnucl == 3)
126 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
128 else
130 fprintf(fplog, "\n");
135 /* TODO The remain call of this function passes non-NULL mass and NULL
136 * mtop, so this routine can be simplified.
138 * The other code path supported doing prediction before the MD loop
139 * started, but even when called, the prediction was always
140 * over-written by a subsequent call in the MD loop, so has been
141 * removed. */
142 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
143 int ns, t_shell s[],
144 real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
146 int i, m, s1, n1, n2, n3;
147 real dt_1, fudge, tm, m1, m2, m3;
148 rvec *ptr;
149 gmx_mtop_atomlookup_t alook = NULL;
150 t_atom *atom;
152 if (mass == NULL)
154 alook = gmx_mtop_atomlookup_init(mtop);
157 /* We introduce a fudge factor for performance reasons: with this choice
158 * the initial force on the shells is about a factor of two lower than
159 * without
161 fudge = 1.0;
163 if (bInit)
165 if (fplog)
167 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
169 ptr = x;
170 dt_1 = 1;
172 else
174 ptr = v;
175 dt_1 = fudge*dt;
178 for (i = 0; (i < ns); i++)
180 s1 = s[i].shell;
181 if (bInit)
183 clear_rvec(x[s1]);
185 switch (s[i].nnucl)
187 case 1:
188 n1 = s[i].nucl1;
189 for (m = 0; (m < DIM); m++)
191 x[s1][m] += ptr[n1][m]*dt_1;
193 break;
194 case 2:
195 n1 = s[i].nucl1;
196 n2 = s[i].nucl2;
197 if (mass)
199 m1 = mass[n1];
200 m2 = mass[n2];
202 else
204 /* Not the correct masses with FE, but it is just a prediction... */
205 gmx_mtop_atomnr_to_atom(alook, n1, &atom);
206 m1 = atom->m;
207 gmx_mtop_atomnr_to_atom(alook, n2, &atom);
208 m2 = atom->m;
210 tm = dt_1/(m1+m2);
211 for (m = 0; (m < DIM); m++)
213 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
215 break;
216 case 3:
217 n1 = s[i].nucl1;
218 n2 = s[i].nucl2;
219 n3 = s[i].nucl3;
220 if (mass)
222 m1 = mass[n1];
223 m2 = mass[n2];
224 m3 = mass[n3];
226 else
228 /* Not the correct masses with FE, but it is just a prediction... */
229 gmx_mtop_atomnr_to_atom(alook, n1, &atom);
230 m1 = atom->m;
231 gmx_mtop_atomnr_to_atom(alook, n2, &atom);
232 m2 = atom->m;
233 gmx_mtop_atomnr_to_atom(alook, n3, &atom);
234 m3 = atom->m;
236 tm = dt_1/(m1+m2+m3);
237 for (m = 0; (m < DIM); m++)
239 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
241 break;
242 default:
243 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
247 if (mass == NULL)
249 gmx_mtop_atomlookup_destroy(alook);
253 /*! \brief Count the different particle types in a system
255 * Routine prints a warning to stderr in case an unknown particle type
256 * is encountered.
257 * \param[in] fplog Print what we have found if not NULL
258 * \param[in] mtop Molecular topology.
259 * \returns Array holding the number of particles of a type
261 static std::array<int, eptNR> countPtypes(FILE *fplog,
262 gmx_mtop_t *mtop)
264 std::array<int, eptNR> nptype = { { 0 } };
265 /* Count number of shells, and find their indices */
266 for (int i = 0; (i < eptNR); i++)
268 nptype[i] = 0;
271 gmx_mtop_atomloop_block_t aloopb = gmx_mtop_atomloop_block_init(mtop);
272 int nmol;
273 t_atom *atom;
274 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
276 switch (atom->ptype)
278 case eptAtom:
279 case eptVSite:
280 case eptShell:
281 nptype[atom->ptype] += nmol;
282 break;
283 default:
284 fprintf(stderr, "Warning unsupported particle type %d in countPtypes",
285 static_cast<int>(atom->ptype));
288 if (fplog)
290 /* Print the number of each particle type */
291 int n = 0;
292 for (const auto &i : nptype)
294 if (i != 0)
296 fprintf(fplog, "There are: %d %ss\n", i, ptype_str[n]);
298 n++;
301 return nptype;
304 gmx_shellfc_t *init_shell_flexcon(FILE *fplog,
305 gmx_mtop_t *mtop, int nflexcon,
306 int nstcalcenergy,
307 bool usingDomainDecomposition)
309 gmx_shellfc_t *shfc;
310 t_shell *shell;
311 int *shell_index = NULL, *at2cg;
312 t_atom *atom;
314 int ns, nshell, nsi;
315 int i, j, type, mb, a_offset, cg, mol, ftype, nra;
316 real qS, alpha;
317 int aS, aN = 0; /* Shell and nucleus */
318 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
319 #define NBT asize(bondtypes)
320 t_iatom *ia;
321 gmx_mtop_atomloop_all_t aloop;
322 gmx_ffparams_t *ffparams;
323 gmx_molblock_t *molb;
324 gmx_moltype_t *molt;
325 t_block *cgs;
327 std::array<int, eptNR> n = countPtypes(fplog, mtop);
328 nshell = n[eptShell];
330 if (nshell == 0 && nflexcon == 0)
332 /* We're not doing shells or flexible constraints */
333 return NULL;
336 snew(shfc, 1);
337 shfc->nflexcon = nflexcon;
339 if (nshell == 0)
341 /* Only flexible constraints, no shells.
342 * Note that make_local_shells() does not need to be called.
344 shfc->nshell = 0;
345 shfc->bPredict = FALSE;
347 return shfc;
350 if (nstcalcenergy != 1)
352 gmx_fatal(FARGS, "You have nstcalcenergy set to a value (%d) that is different from 1.\nThis is not supported in combination with shell particles.\nPlease make a new tpr file.", nstcalcenergy);
354 if (usingDomainDecomposition)
356 gmx_fatal(FARGS, "Shell particles are not implemented with domain decomposition, use a single rank");
359 /* We have shells: fill the shell data structure */
361 /* Global system sized array, this should be avoided */
362 snew(shell_index, mtop->natoms);
364 aloop = gmx_mtop_atomloop_all_init(mtop);
365 nshell = 0;
366 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
368 if (atom->ptype == eptShell)
370 shell_index[i] = nshell++;
374 snew(shell, nshell);
376 /* Initiate the shell structures */
377 for (i = 0; (i < nshell); i++)
379 shell[i].shell = -1;
380 shell[i].nnucl = 0;
381 shell[i].nucl1 = -1;
382 shell[i].nucl2 = -1;
383 shell[i].nucl3 = -1;
384 /* shell[i].bInterCG=FALSE; */
385 shell[i].k_1 = 0;
386 shell[i].k = 0;
389 ffparams = &mtop->ffparams;
391 /* Now fill the structures */
392 shfc->bInterCG = FALSE;
393 ns = 0;
394 a_offset = 0;
395 for (mb = 0; mb < mtop->nmolblock; mb++)
397 molb = &mtop->molblock[mb];
398 molt = &mtop->moltype[molb->type];
400 cgs = &molt->cgs;
401 snew(at2cg, molt->atoms.nr);
402 for (cg = 0; cg < cgs->nr; cg++)
404 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
406 at2cg[i] = cg;
410 atom = molt->atoms.atom;
411 for (mol = 0; mol < molb->nmol; mol++)
413 for (j = 0; (j < NBT); j++)
415 ia = molt->ilist[bondtypes[j]].iatoms;
416 for (i = 0; (i < molt->ilist[bondtypes[j]].nr); )
418 type = ia[0];
419 ftype = ffparams->functype[type];
420 nra = interaction_function[ftype].nratoms;
422 /* Check whether we have a bond with a shell */
423 aS = -1;
425 switch (bondtypes[j])
427 case F_BONDS:
428 case F_HARMONIC:
429 case F_CUBICBONDS:
430 case F_POLARIZATION:
431 case F_ANHARM_POL:
432 if (atom[ia[1]].ptype == eptShell)
434 aS = ia[1];
435 aN = ia[2];
437 else if (atom[ia[2]].ptype == eptShell)
439 aS = ia[2];
440 aN = ia[1];
442 break;
443 case F_WATER_POL:
444 aN = ia[4]; /* Dummy */
445 aS = ia[5]; /* Shell */
446 break;
447 default:
448 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
451 if (aS != -1)
453 qS = atom[aS].q;
455 /* Check whether one of the particles is a shell... */
456 nsi = shell_index[a_offset+aS];
457 if ((nsi < 0) || (nsi >= nshell))
459 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
460 nsi, nshell, aS);
462 if (shell[nsi].shell == -1)
464 shell[nsi].shell = a_offset + aS;
465 ns++;
467 else if (shell[nsi].shell != a_offset+aS)
469 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
472 if (shell[nsi].nucl1 == -1)
474 shell[nsi].nucl1 = a_offset + aN;
476 else if (shell[nsi].nucl2 == -1)
478 shell[nsi].nucl2 = a_offset + aN;
480 else if (shell[nsi].nucl3 == -1)
482 shell[nsi].nucl3 = a_offset + aN;
484 else
486 if (fplog)
488 pr_shell(fplog, ns, shell);
490 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
492 if (at2cg[aS] != at2cg[aN])
494 /* shell[nsi].bInterCG = TRUE; */
495 shfc->bInterCG = TRUE;
498 switch (bondtypes[j])
500 case F_BONDS:
501 case F_HARMONIC:
502 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
503 break;
504 case F_CUBICBONDS:
505 shell[nsi].k += ffparams->iparams[type].cubic.kb;
506 break;
507 case F_POLARIZATION:
508 case F_ANHARM_POL:
509 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
511 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);
513 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/
514 ffparams->iparams[type].polarize.alpha;
515 break;
516 case F_WATER_POL:
517 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
519 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);
521 alpha = (ffparams->iparams[type].wpol.al_x+
522 ffparams->iparams[type].wpol.al_y+
523 ffparams->iparams[type].wpol.al_z)/3.0;
524 shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/alpha;
525 break;
526 default:
527 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
529 shell[nsi].nnucl++;
531 ia += nra+1;
532 i += nra+1;
535 a_offset += molt->atoms.nr;
537 /* Done with this molecule type */
538 sfree(at2cg);
541 /* Verify whether it's all correct */
542 if (ns != nshell)
544 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
547 for (i = 0; (i < ns); i++)
549 shell[i].k_1 = 1.0/shell[i].k;
552 if (debug)
554 pr_shell(debug, ns, shell);
558 shfc->nshell_gl = ns;
559 shfc->shell_gl = shell;
560 shfc->shell_index_gl = shell_index;
562 shfc->bPredict = (getenv("GMX_NOPREDICT") == NULL);
563 shfc->bRequireInit = FALSE;
564 if (!shfc->bPredict)
566 if (fplog)
568 fprintf(fplog, "\nWill never predict shell positions\n");
571 else
573 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != NULL);
574 if (shfc->bRequireInit && fplog)
576 fprintf(fplog, "\nWill always initiate shell positions\n");
580 if (shfc->bPredict)
582 if (shfc->bInterCG)
584 if (fplog)
586 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");
588 /* Prediction improves performance, so we should implement either:
589 * 1. communication for the atoms needed for prediction
590 * 2. prediction using the velocities of shells; currently the
591 * shell velocities are zeroed, it's a bit tricky to keep
592 * track of the shell displacements and thus the velocity.
594 shfc->bPredict = FALSE;
598 return shfc;
601 void make_local_shells(t_commrec *cr, t_mdatoms *md,
602 gmx_shellfc_t *shfc)
604 t_shell *shell;
605 int a0, a1, *ind, nshell, i;
606 gmx_domdec_t *dd = NULL;
608 if (DOMAINDECOMP(cr))
610 dd = cr->dd;
611 a0 = 0;
612 a1 = dd->nat_home;
614 else
616 /* Single node: we need all shells, just copy the pointer */
617 shfc->nshell = shfc->nshell_gl;
618 shfc->shell = shfc->shell_gl;
620 return;
623 ind = shfc->shell_index_gl;
625 nshell = 0;
626 shell = shfc->shell;
627 for (i = a0; i < a1; i++)
629 if (md->ptype[i] == eptShell)
631 if (nshell+1 > shfc->shell_nalloc)
633 shfc->shell_nalloc = over_alloc_dd(nshell+1);
634 srenew(shell, shfc->shell_nalloc);
636 if (dd)
638 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
640 else
642 shell[nshell] = shfc->shell_gl[ind[i]];
645 /* With inter-cg shells we can no do shell prediction,
646 * so we do not need the nuclei numbers.
648 if (!shfc->bInterCG)
650 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
651 if (shell[nshell].nnucl > 1)
653 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
655 if (shell[nshell].nnucl > 2)
657 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
660 shell[nshell].shell = i;
661 nshell++;
665 shfc->nshell = nshell;
666 shfc->shell = shell;
669 static void do_1pos(rvec xnew, rvec xold, rvec f, real step)
671 real xo, yo, zo;
672 real dx, dy, dz;
674 xo = xold[XX];
675 yo = xold[YY];
676 zo = xold[ZZ];
678 dx = f[XX]*step;
679 dy = f[YY]*step;
680 dz = f[ZZ]*step;
682 xnew[XX] = xo+dx;
683 xnew[YY] = yo+dy;
684 xnew[ZZ] = zo+dz;
687 static void do_1pos3(rvec xnew, rvec xold, rvec f, rvec step)
689 real xo, yo, zo;
690 real dx, dy, dz;
692 xo = xold[XX];
693 yo = xold[YY];
694 zo = xold[ZZ];
696 dx = f[XX]*step[XX];
697 dy = f[YY]*step[YY];
698 dz = f[ZZ]*step[ZZ];
700 xnew[XX] = xo+dx;
701 xnew[YY] = yo+dy;
702 xnew[ZZ] = zo+dz;
705 static void directional_sd(rvec xold[], rvec xnew[], rvec acc_dir[],
706 int start, int homenr, real step)
708 int i;
710 for (i = start; i < homenr; i++)
712 do_1pos(xnew[i], xold[i], acc_dir[i], step);
716 static void shell_pos_sd(rvec xcur[], rvec xnew[], rvec f[],
717 int ns, t_shell s[], int count)
719 const real step_scale_min = 0.8,
720 step_scale_increment = 0.2,
721 step_scale_max = 1.2,
722 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
723 int i, shell, d;
724 real dx, df, k_est;
725 const real zero = 0;
726 #ifdef PRINT_STEP
727 real step_min, step_max;
729 step_min = 1e30;
730 step_max = 0;
731 #endif
732 for (i = 0; (i < ns); i++)
734 shell = s[i].shell;
735 if (count == 1)
737 for (d = 0; d < DIM; d++)
739 s[i].step[d] = s[i].k_1;
740 #ifdef PRINT_STEP
741 step_min = std::min(step_min, s[i].step[d]);
742 step_max = std::max(step_max, s[i].step[d]);
743 #endif
746 else
748 for (d = 0; d < DIM; d++)
750 dx = xcur[shell][d] - s[i].xold[d];
751 df = f[shell][d] - s[i].fold[d];
752 /* -dx/df gets used to generate an interpolated value, but would
753 * cause a NaN if df were binary-equal to zero. Values close to
754 * zero won't cause problems (because of the min() and max()), so
755 * just testing for binary inequality is OK. */
756 if (zero != df)
758 k_est = -dx/df;
759 /* Scale the step size by a factor interpolated from
760 * step_scale_min to step_scale_max, as k_est goes from 0 to
761 * step_scale_multiple * s[i].step[d] */
762 s[i].step[d] =
763 step_scale_min * s[i].step[d] +
764 step_scale_increment * std::min(step_scale_multiple * s[i].step[d], std::max(k_est, zero));
766 else
768 /* Here 0 == df */
769 if (gmx_numzero(dx)) /* 0 == dx */
771 /* Likely this will never happen, but if it does just
772 * don't scale the step. */
774 else /* 0 != dx */
776 s[i].step[d] *= step_scale_max;
779 #ifdef PRINT_STEP
780 step_min = std::min(step_min, s[i].step[d]);
781 step_max = std::max(step_max, s[i].step[d]);
782 #endif
785 copy_rvec(xcur[shell], s[i].xold);
786 copy_rvec(f[shell], s[i].fold);
788 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
790 if (gmx_debug_at)
792 fprintf(debug, "shell[%d] = %d\n", i, shell);
793 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
794 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
795 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
796 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
799 #ifdef PRINT_STEP
800 printf("step %.3e %.3e\n", step_min, step_max);
801 #endif
804 static void decrease_step_size(int nshell, t_shell s[])
806 int i;
808 for (i = 0; i < nshell; i++)
810 svmul(0.8, s[i].step, s[i].step);
814 static void print_epot(FILE *fp, gmx_int64_t mdstep, int count, real epot, real df,
815 int ndir, real sf_dir)
817 char buf[22];
819 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
820 gmx_step_str(mdstep, buf), count, epot, df);
821 if (ndir)
823 fprintf(fp, ", dir. rmsF: %6.2e\n", std::sqrt(sf_dir/ndir));
825 else
827 fprintf(fp, "\n");
832 static real rms_force(t_commrec *cr, rvec f[], int ns, t_shell s[],
833 int ndir, real *sf_dir, real *Epot)
835 int i, shell, ntot;
836 double buf[4];
838 buf[0] = *sf_dir;
839 for (i = 0; i < ns; i++)
841 shell = s[i].shell;
842 buf[0] += norm2(f[shell]);
844 ntot = ns;
846 if (PAR(cr))
848 buf[1] = ntot;
849 buf[2] = *sf_dir;
850 buf[3] = *Epot;
851 gmx_sumd(4, buf, cr);
852 ntot = (int)(buf[1] + 0.5);
853 *sf_dir = buf[2];
854 *Epot = buf[3];
856 ntot += ndir;
858 return (ntot ? std::sqrt(buf[0]/ntot) : 0);
861 static void check_pbc(FILE *fp, rvec x[], int shell)
863 int m, now;
865 now = shell-4;
866 for (m = 0; (m < DIM); m++)
868 if (fabs(x[shell][m]-x[now][m]) > 0.3)
870 pr_rvecs(fp, 0, "SHELL-X", x+now, 5);
871 break;
876 static void dump_shells(FILE *fp, rvec x[], rvec f[], real ftol, int ns, t_shell s[])
878 int i, shell;
879 real ft2, ff2;
881 ft2 = gmx::square(ftol);
883 for (i = 0; (i < ns); i++)
885 shell = s[i].shell;
886 ff2 = iprod(f[shell], f[shell]);
887 if (ff2 > ft2)
889 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
890 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], std::sqrt(ff2));
892 check_pbc(fp, x, shell);
896 static void init_adir(FILE *log, gmx_shellfc_t *shfc,
897 gmx_constr_t constr, t_idef *idef, t_inputrec *ir,
898 t_commrec *cr, int dd_ac1,
899 gmx_int64_t step, t_mdatoms *md, int start, int end,
900 rvec *x_old, rvec *x_init, rvec *x,
901 rvec *f, rvec *acc_dir,
902 gmx_bool bMolPBC, matrix box,
903 real *lambda, real *dvdlambda, t_nrnb *nrnb)
905 rvec *xnold, *xnew;
906 double dt, w_dt;
907 int n, d;
908 unsigned short *ptype;
910 if (DOMAINDECOMP(cr))
912 n = dd_ac1;
914 else
916 n = end - start;
918 if (n > shfc->adir_nalloc)
920 shfc->adir_nalloc = over_alloc_dd(n);
921 srenew(shfc->adir_xnold, shfc->adir_nalloc);
922 srenew(shfc->adir_xnew, shfc->adir_nalloc);
924 xnold = shfc->adir_xnold;
925 xnew = shfc->adir_xnew;
927 ptype = md->ptype;
929 dt = ir->delta_t;
931 /* Does NOT work with freeze or acceleration groups (yet) */
932 for (n = start; n < end; n++)
934 w_dt = md->invmass[n]*dt;
936 for (d = 0; d < DIM; d++)
938 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
940 xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
941 xnew[n-start][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
943 else
945 xnold[n-start][d] = x[n][d];
946 xnew[n-start][d] = x[n][d];
950 constrain(log, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
951 x, xnold-start, NULL, bMolPBC, box,
952 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
953 NULL, NULL, nrnb, econqCoord);
954 constrain(log, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
955 x, xnew-start, NULL, bMolPBC, box,
956 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
957 NULL, NULL, nrnb, econqCoord);
959 for (n = start; n < end; n++)
961 for (d = 0; d < DIM; d++)
963 xnew[n-start][d] =
964 -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/gmx::square(dt)
965 - f[n][d]*md->invmass[n];
967 clear_rvec(acc_dir[n]);
970 /* Project the acceleration on the old bond directions */
971 constrain(log, FALSE, FALSE, constr, idef, ir, cr, step, 0, 1.0, md,
972 x_old, xnew-start, acc_dir, bMolPBC, box,
973 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
974 NULL, NULL, nrnb, econqDeriv_FlexCon);
977 void relax_shell_flexcon(FILE *fplog, t_commrec *cr, gmx_bool bVerbose,
978 gmx_int64_t mdstep, t_inputrec *inputrec,
979 gmx_bool bDoNS, int force_flags,
980 gmx_localtop_t *top,
981 gmx_constr_t constr,
982 gmx_enerdata_t *enerd, t_fcdata *fcd,
983 t_state *state, rvec f[],
984 tensor force_vir,
985 t_mdatoms *md,
986 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
987 t_graph *graph,
988 gmx_groups_t *groups,
989 gmx_shellfc_t *shfc,
990 t_forcerec *fr,
991 gmx_bool bBornRadii,
992 double t, rvec mu_tot,
993 gmx_vsite_t *vsite,
994 FILE *fp_field)
996 int nshell;
997 t_shell *shell;
998 t_idef *idef;
999 rvec *pos[2], *force[2], *acc_dir = NULL, *x_old = NULL;
1000 real Epot[2], df[2];
1001 real sf_dir, invdt;
1002 real ftol, dum = 0;
1003 char sbuf[22];
1004 gmx_bool bCont, bInit, bConverged;
1005 int nat, dd_ac0, dd_ac1 = 0, i;
1006 int start = 0, homenr = md->homenr, end = start+homenr, cg0, cg1;
1007 int nflexcon, number_steps, d, Min = 0, count = 0;
1008 #define Try (1-Min) /* At start Try = 1 */
1010 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
1011 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
1012 ftol = inputrec->em_tol;
1013 number_steps = inputrec->niter;
1014 nshell = shfc->nshell;
1015 shell = shfc->shell;
1016 nflexcon = shfc->nflexcon;
1018 idef = &top->idef;
1020 if (DOMAINDECOMP(cr))
1022 nat = dd_natoms_vsite(cr->dd);
1023 if (nflexcon > 0)
1025 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
1026 nat = std::max(nat, dd_ac1);
1029 else
1031 nat = state->natoms;
1034 if (nat > shfc->x_nalloc)
1036 /* Allocate local arrays */
1037 shfc->x_nalloc = over_alloc_dd(nat);
1038 for (i = 0; (i < 2); i++)
1040 srenew(shfc->x[i], shfc->x_nalloc);
1041 srenew(shfc->f[i], shfc->x_nalloc);
1044 for (i = 0; (i < 2); i++)
1046 pos[i] = shfc->x[i];
1047 force[i] = shfc->f[i];
1050 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
1052 /* This is the only time where the coordinates are used
1053 * before do_force is called, which normally puts all
1054 * charge groups in the box.
1056 if (inputrec->cutoff_scheme == ecutsVERLET)
1058 put_atoms_in_box_omp(fr->ePBC, state->box, md->homenr, state->x);
1060 else
1062 cg0 = 0;
1063 cg1 = top->cgs.nr;
1064 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1065 &(top->cgs), state->x, fr->cg_cm);
1068 if (graph)
1070 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
1074 /* After this all coordinate arrays will contain whole charge groups */
1075 if (graph)
1077 shift_self(graph, state->box, state->x);
1080 if (nflexcon)
1082 if (nat > shfc->flex_nalloc)
1084 shfc->flex_nalloc = over_alloc_dd(nat);
1085 srenew(shfc->acc_dir, shfc->flex_nalloc);
1086 srenew(shfc->x_old, shfc->flex_nalloc);
1088 acc_dir = shfc->acc_dir;
1089 x_old = shfc->x_old;
1090 for (i = 0; i < homenr; i++)
1092 for (d = 0; d < DIM; d++)
1094 shfc->x_old[i][d] =
1095 state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
1100 /* Do a prediction of the shell positions */
1101 if (shfc->bPredict && !bCont)
1103 predict_shells(fplog, state->x, state->v, inputrec->delta_t, nshell, shell,
1104 md->massT, NULL, bInit);
1107 /* do_force expected the charge groups to be in the box */
1108 if (graph)
1110 unshift_self(graph, state->box, state->x);
1113 /* Calculate the forces first time around */
1114 if (gmx_debug_at)
1116 pr_rvecs(debug, 0, "x b4 do_force", state->x + start, homenr);
1118 do_force(fplog, cr, inputrec, mdstep, nrnb, wcycle, top, groups,
1119 state->box, state->x, &state->hist,
1120 force[Min], force_vir, md, enerd, fcd,
1121 state->lambda, graph,
1122 fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1123 (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
1125 sf_dir = 0;
1126 if (nflexcon)
1128 init_adir(fplog, shfc,
1129 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1130 shfc->x_old-start, state->x, state->x, force[Min],
1131 shfc->acc_dir-start,
1132 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1134 for (i = start; i < end; i++)
1136 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
1140 Epot[Min] = enerd->term[F_EPOT];
1142 df[Min] = rms_force(cr, shfc->f[Min], nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1143 df[Try] = 0;
1144 if (debug)
1146 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1149 if (gmx_debug_at)
1151 pr_rvecs(debug, 0, "force0", force[Min], md->nr);
1154 if (nshell+nflexcon > 0)
1156 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1157 * shell positions are updated, therefore the other particles must
1158 * be set here.
1160 memcpy(pos[Min], state->x, nat*sizeof(state->x[0]));
1161 memcpy(pos[Try], state->x, nat*sizeof(state->x[0]));
1164 if (bVerbose && MASTER(cr))
1166 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1169 if (debug)
1171 fprintf(debug, "%17s: %14.10e\n",
1172 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1173 fprintf(debug, "%17s: %14.10e\n",
1174 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1175 fprintf(debug, "%17s: %14.10e\n",
1176 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1177 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1180 /* First check whether we should do shells, or whether the force is
1181 * low enough even without minimization.
1183 bConverged = (df[Min] < ftol);
1185 for (count = 1; (!(bConverged) && (count < number_steps)); count++)
1187 if (vsite)
1189 construct_vsites(vsite, pos[Min], inputrec->delta_t, state->v,
1190 idef->iparams, idef->il,
1191 fr->ePBC, fr->bMolPBC, cr, state->box);
1194 if (nflexcon)
1196 init_adir(fplog, shfc,
1197 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1198 x_old-start, state->x, pos[Min], force[Min], acc_dir-start,
1199 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1201 directional_sd(pos[Min], pos[Try], acc_dir-start, start, end,
1202 fr->fc_stepsize);
1205 /* New positions, Steepest descent */
1206 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1208 /* do_force expected the charge groups to be in the box */
1209 if (graph)
1211 unshift_self(graph, state->box, pos[Try]);
1214 if (gmx_debug_at)
1216 pr_rvecs(debug, 0, "RELAX: pos[Min] ", pos[Min] + start, homenr);
1217 pr_rvecs(debug, 0, "RELAX: pos[Try] ", pos[Try] + start, homenr);
1219 /* Try the new positions */
1220 do_force(fplog, cr, inputrec, 1, nrnb, wcycle,
1221 top, groups, state->box, pos[Try], &state->hist,
1222 force[Try], force_vir,
1223 md, enerd, fcd, state->lambda, graph,
1224 fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1225 force_flags);
1227 if (gmx_debug_at)
1229 pr_rvecs(debug, 0, "RELAX: force[Min]", force[Min] + start, homenr);
1230 pr_rvecs(debug, 0, "RELAX: force[Try]", force[Try] + start, homenr);
1232 sf_dir = 0;
1233 if (nflexcon)
1235 init_adir(fplog, shfc,
1236 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1237 x_old-start, state->x, pos[Try], force[Try], acc_dir-start,
1238 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1240 for (i = start; i < end; i++)
1242 sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
1246 Epot[Try] = enerd->term[F_EPOT];
1248 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1250 if (debug)
1252 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1255 if (debug)
1257 if (gmx_debug_at)
1259 pr_rvecs(debug, 0, "F na do_force", force[Try] + start, homenr);
1261 if (gmx_debug_at)
1263 fprintf(debug, "SHELL ITER %d\n", count);
1264 dump_shells(debug, pos[Try], force[Try], ftol, nshell, shell);
1268 if (bVerbose && MASTER(cr))
1270 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1273 bConverged = (df[Try] < ftol);
1275 if ((df[Try] < df[Min]))
1277 if (debug)
1279 fprintf(debug, "Swapping Min and Try\n");
1281 if (nflexcon)
1283 /* Correct the velocities for the flexible constraints */
1284 invdt = 1/inputrec->delta_t;
1285 for (i = start; i < end; i++)
1287 for (d = 0; d < DIM; d++)
1289 state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1293 Min = Try;
1295 else
1297 decrease_step_size(nshell, shell);
1300 shfc->numForceEvaluations += count;
1301 if (bConverged)
1303 shfc->numConvergedIterations++;
1305 if (MASTER(cr) && !(bConverged))
1307 /* Note that the energies and virial are incorrect when not converged */
1308 if (fplog)
1310 fprintf(fplog,
1311 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1312 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1314 fprintf(stderr,
1315 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1316 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1319 /* Copy back the coordinates and the forces */
1320 memcpy(state->x, pos[Min], nat*sizeof(state->x[0]));
1321 memcpy(f, force[Min], nat*sizeof(f[0]));
1324 void done_shellfc(FILE *fplog, gmx_shellfc_t *shfc, gmx_int64_t numSteps)
1326 if (shfc && fplog && numSteps > 0)
1328 double numStepsAsDouble = static_cast<double>(numSteps);
1329 fprintf(fplog, "Fraction of iterations that converged: %.2f %%\n",
1330 (shfc->numConvergedIterations*100.0)/numStepsAsDouble);
1331 fprintf(fplog, "Average number of force evaluations per MD step: %.2f\n\n",
1332 shfc->numForceEvaluations/numStepsAsDouble);
1335 // TODO Deallocate memory in shfc