3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2008, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
32 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
41 #include "gmx_fatal.h"
56 #include "mtop_util.h"
57 #include "chargegroup.h"
63 atom_id shell; /* The shell id */
64 atom_id nucl1, nucl2, nucl3; /* The nuclei connected to the shell */
65 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
66 real k; /* force constant */
67 real k_1; /* 1 over force constant */
73 typedef struct gmx_shellfc {
74 int nshell_gl; /* The number of shells in the system */
75 t_shell *shell_gl; /* All the shells (for DD only) */
76 int *shell_index_gl; /* Global shell index (for DD only) */
77 gmx_bool bInterCG; /* Are there inter charge-group shells? */
78 int nshell; /* The number of local shells */
79 t_shell *shell; /* The local shells */
80 int shell_nalloc; /* The allocation size of shell */
81 gmx_bool bPredict; /* Predict shell positions */
82 gmx_bool bRequireInit; /* Require initialization of shell positions */
83 int nflexcon; /* The number of flexible constraints */
84 rvec *x[2]; /* Array for iterative minimization */
85 rvec *f[2]; /* Array for iterative minimization */
86 int x_nalloc; /* The allocation size of x and f */
87 rvec *acc_dir; /* Acceleration direction for flexcon */
88 rvec *x_old; /* Old coordinates for flexcon */
89 int flex_nalloc; /* The allocation size of acc_dir and x_old */
90 rvec *adir_xnold; /* Work space for init_adir */
91 rvec *adir_xnew; /* Work space for init_adir */
92 int adir_nalloc; /* Work space for init_adir */
96 static void pr_shell(FILE *fplog, int ns, t_shell s[])
100 fprintf(fplog, "SHELL DATA\n");
101 fprintf(fplog, "%5s %8s %5s %5s %5s\n",
102 "Shell", "Force k", "Nucl1", "Nucl2", "Nucl3");
103 for (i = 0; (i < ns); i++)
105 fprintf(fplog, "%5d %8.3f %5d", s[i].shell, 1.0/s[i].k_1, s[i].nucl1);
108 fprintf(fplog, " %5d\n", s[i].nucl2);
110 else if (s[i].nnucl == 3)
112 fprintf(fplog, " %5d %5d\n", s[i].nucl2, s[i].nucl3);
116 fprintf(fplog, "\n");
121 static void predict_shells(FILE *fplog, rvec x[], rvec v[], real dt,
123 real mass[], gmx_mtop_t *mtop, gmx_bool bInit)
125 int i, m, s1, n1, n2, n3;
126 real dt_1, dt_2, dt_3, fudge, tm, m1, m2, m3;
128 gmx_mtop_atomlookup_t alook = NULL;
133 alook = gmx_mtop_atomlookup_init(mtop);
136 /* We introduce a fudge factor for performance reasons: with this choice
137 * the initial force on the shells is about a factor of two lower than
146 fprintf(fplog, "RELAX: Using prediction for initial shell placement\n");
157 for (i = 0; (i < ns); i++)
168 for (m = 0; (m < DIM); m++)
170 x[s1][m] += ptr[n1][m]*dt_1;
183 /* Not the correct masses with FE, but it is just a prediction... */
188 for (m = 0; (m < DIM); m++)
190 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
205 /* Not the correct masses with FE, but it is just a prediction... */
206 gmx_mtop_atomnr_to_atom(alook, n1, &atom);
208 gmx_mtop_atomnr_to_atom(alook, n2, &atom);
210 gmx_mtop_atomnr_to_atom(alook, n3, &atom);
213 tm = dt_1/(m1+m2+m3);
214 for (m = 0; (m < DIM); m++)
216 x[s1][m] += (m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
220 gmx_fatal(FARGS, "Shell %d has %d nuclei!", i, s[i].nnucl);
226 gmx_mtop_atomlookup_destroy(alook);
230 gmx_shellfc_t init_shell_flexcon(FILE *fplog,
231 gmx_mtop_t *mtop, int nflexcon,
234 struct gmx_shellfc *shfc;
236 int *shell_index = NULL, *at2cg;
238 int n[eptNR], ns, nshell, nsi;
239 int i, j, nmol, type, mb, mt, a_offset, cg, mol, ftype, nra;
241 int aS, aN = 0; /* Shell and nucleus */
242 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
243 #define NBT asize(bondtypes)
245 gmx_mtop_atomloop_block_t aloopb;
246 gmx_mtop_atomloop_all_t aloop;
247 gmx_ffparams_t *ffparams;
248 gmx_molblock_t *molb;
252 /* Count number of shells, and find their indices */
253 for (i = 0; (i < eptNR); i++)
258 aloopb = gmx_mtop_atomloop_block_init(mtop);
259 while (gmx_mtop_atomloop_block_next(aloopb, &atom, &nmol))
261 n[atom->ptype] += nmol;
266 /* Print the number of each particle type */
267 for (i = 0; (i < eptNR); i++)
271 fprintf(fplog, "There are: %d %ss\n", n[i], ptype_str[i]);
276 nshell = n[eptShell];
278 if (nshell == 0 && nflexcon == 0)
284 shfc->nflexcon = nflexcon;
291 /* We have shells: fill the shell data structure */
293 /* Global system sized array, this should be avoided */
294 snew(shell_index, mtop->natoms);
296 aloop = gmx_mtop_atomloop_all_init(mtop);
298 while (gmx_mtop_atomloop_all_next(aloop, &i, &atom))
300 if (atom->ptype == eptShell)
302 shell_index[i] = nshell++;
308 /* Initiate the shell structures */
309 for (i = 0; (i < nshell); i++)
311 shell[i].shell = NO_ATID;
313 shell[i].nucl1 = NO_ATID;
314 shell[i].nucl2 = NO_ATID;
315 shell[i].nucl3 = NO_ATID;
316 /* shell[i].bInterCG=FALSE; */
321 ffparams = &mtop->ffparams;
323 /* Now fill the structures */
324 shfc->bInterCG = FALSE;
327 for (mb = 0; mb < mtop->nmolblock; mb++)
329 molb = &mtop->molblock[mb];
330 molt = &mtop->moltype[molb->type];
333 snew(at2cg, molt->atoms.nr);
334 for (cg = 0; cg < cgs->nr; cg++)
336 for (i = cgs->index[cg]; i < cgs->index[cg+1]; i++)
342 atom = molt->atoms.atom;
343 for (mol = 0; mol < molb->nmol; mol++)
345 for (j = 0; (j < NBT); j++)
347 ia = molt->ilist[bondtypes[j]].iatoms;
348 for (i = 0; (i < molt->ilist[bondtypes[j]].nr); )
351 ftype = ffparams->functype[type];
352 nra = interaction_function[ftype].nratoms;
354 /* Check whether we have a bond with a shell */
357 switch (bondtypes[j])
364 if (atom[ia[1]].ptype == eptShell)
369 else if (atom[ia[2]].ptype == eptShell)
376 aN = ia[4]; /* Dummy */
377 aS = ia[5]; /* Shell */
380 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
387 /* Check whether one of the particles is a shell... */
388 nsi = shell_index[a_offset+aS];
389 if ((nsi < 0) || (nsi >= nshell))
391 gmx_fatal(FARGS, "nsi is %d should be within 0 - %d. aS = %d",
394 if (shell[nsi].shell == NO_ATID)
396 shell[nsi].shell = a_offset + aS;
399 else if (shell[nsi].shell != a_offset+aS)
401 gmx_fatal(FARGS, "Weird stuff in %s, %d", __FILE__, __LINE__);
404 if (shell[nsi].nucl1 == NO_ATID)
406 shell[nsi].nucl1 = a_offset + aN;
408 else if (shell[nsi].nucl2 == NO_ATID)
410 shell[nsi].nucl2 = a_offset + aN;
412 else if (shell[nsi].nucl3 == NO_ATID)
414 shell[nsi].nucl3 = a_offset + aN;
420 pr_shell(fplog, ns, shell);
422 gmx_fatal(FARGS, "Can not handle more than three bonds per shell\n");
424 if (at2cg[aS] != at2cg[aN])
426 /* shell[nsi].bInterCG = TRUE; */
427 shfc->bInterCG = TRUE;
430 switch (bondtypes[j])
434 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
437 shell[nsi].k += ffparams->iparams[type].cubic.kb;
441 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
443 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 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/
446 ffparams->iparams[type].polarize.alpha;
449 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
451 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 alpha = (ffparams->iparams[type].wpol.al_x+
454 ffparams->iparams[type].wpol.al_y+
455 ffparams->iparams[type].wpol.al_z)/3.0;
456 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/alpha;
459 gmx_fatal(FARGS, "Death Horror: %s, %d", __FILE__, __LINE__);
467 a_offset += molt->atoms.nr;
469 /* Done with this molecule type */
473 /* Verify whether it's all correct */
476 gmx_fatal(FARGS, "Something weird with shells. They may not be bonded to something");
479 for (i = 0; (i < ns); i++)
481 shell[i].k_1 = 1.0/shell[i].k;
486 pr_shell(debug, ns, shell);
490 shfc->nshell_gl = ns;
491 shfc->shell_gl = shell;
492 shfc->shell_index_gl = shell_index;
494 shfc->bPredict = (getenv("GMX_NOPREDICT") == NULL);
495 shfc->bRequireInit = FALSE;
500 fprintf(fplog, "\nWill never predict shell positions\n");
505 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != NULL);
506 if (shfc->bRequireInit && fplog)
508 fprintf(fplog, "\nWill always initiate shell positions\n");
516 predict_shells(fplog, x, NULL, 0, shfc->nshell_gl, shfc->shell_gl,
524 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 shfc->bPredict = FALSE;
533 void make_local_shells(t_commrec *cr, t_mdatoms *md,
534 struct gmx_shellfc *shfc)
537 int a0, a1, *ind, nshell, i;
538 gmx_domdec_t *dd = NULL;
542 if (DOMAINDECOMP(cr))
550 pd_at_range(cr, &a0, &a1);
555 /* Single node: we need all shells, just copy the pointer */
556 shfc->nshell = shfc->nshell_gl;
557 shfc->shell = shfc->shell_gl;
562 ind = shfc->shell_index_gl;
566 for (i = a0; i < a1; i++)
568 if (md->ptype[i] == eptShell)
570 if (nshell+1 > shfc->shell_nalloc)
572 shfc->shell_nalloc = over_alloc_dd(nshell+1);
573 srenew(shell, shfc->shell_nalloc);
577 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
581 shell[nshell] = shfc->shell_gl[ind[i]];
583 /* With inter-cg shells we can no do shell prediction,
584 * so we do not need the nuclei numbers.
588 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
589 if (shell[nshell].nnucl > 1)
591 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
593 if (shell[nshell].nnucl > 2)
595 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
598 shell[nshell].shell = i;
603 shfc->nshell = nshell;
607 static void do_1pos(rvec xnew, rvec xold, rvec f, real step)
625 static void do_1pos3(rvec xnew, rvec xold, rvec f, rvec step)
643 static void directional_sd(rvec xold[], rvec xnew[], rvec acc_dir[],
644 int start, int homenr, real step)
648 for (i = start; i < homenr; i++)
650 do_1pos(xnew[i], xold[i], acc_dir[i], step);
654 static void shell_pos_sd(rvec xcur[], rvec xnew[], rvec f[],
655 int ns, t_shell s[], int count)
657 const real step_scale_min = 0.8,
658 step_scale_increment = 0.2,
659 step_scale_max = 1.2,
660 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
664 real step_min, step_max;
669 for (i = 0; (i < ns); i++)
674 for (d = 0; d < DIM; d++)
676 s[i].step[d] = s[i].k_1;
678 step_min = min(step_min, s[i].step[d]);
679 step_max = max(step_max, s[i].step[d]);
685 for (d = 0; d < DIM; d++)
687 dx = xcur[shell][d] - s[i].xold[d];
688 df = f[shell][d] - s[i].fold[d];
689 /* -dx/df gets used to generate an interpolated value, but would
690 * cause a NaN if df were binary-equal to zero. Values close to
691 * zero won't cause problems (because of the min() and max()), so
692 * just testing for binary inequality is OK. */
696 /* Scale the step size by a factor interpolated from
697 * step_scale_min to step_scale_max, as k_est goes from 0 to
698 * step_scale_multiple * s[i].step[d] */
700 step_scale_min * s[i].step[d] +
701 step_scale_increment * min(step_scale_multiple * s[i].step[d], max(k_est, 0));
706 if (gmx_numzero(dx)) /* 0 == dx */
708 /* Likely this will never happen, but if it does just
709 * don't scale the step. */
713 s[i].step[d] *= step_scale_max;
717 step_min = min(step_min, s[i].step[d]);
718 step_max = max(step_max, s[i].step[d]);
722 copy_rvec(xcur[shell], s[i].xold);
723 copy_rvec(f[shell], s[i].fold);
725 do_1pos3(xnew[shell], xcur[shell], f[shell], s[i].step);
729 fprintf(debug, "shell[%d] = %d\n", i, shell);
730 pr_rvec(debug, 0, "fshell", f[shell], DIM, TRUE);
731 pr_rvec(debug, 0, "xold", xcur[shell], DIM, TRUE);
732 pr_rvec(debug, 0, "step", s[i].step, DIM, TRUE);
733 pr_rvec(debug, 0, "xnew", xnew[shell], DIM, TRUE);
737 printf("step %.3e %.3e\n", step_min, step_max);
741 static void decrease_step_size(int nshell, t_shell s[])
745 for (i = 0; i < nshell; i++)
747 svmul(0.8, s[i].step, s[i].step);
751 static void print_epot(FILE *fp, gmx_large_int_t mdstep, int count, real epot, real df,
752 int ndir, real sf_dir)
756 fprintf(fp, "MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
757 gmx_step_str(mdstep, buf), count, epot, df);
760 fprintf(fp, ", dir. rmsF: %6.2e\n", sqrt(sf_dir/ndir));
769 static real rms_force(t_commrec *cr, rvec f[], int ns, t_shell s[],
770 int ndir, real *sf_dir, real *Epot)
776 for (i = 0; i < ns; i++)
779 buf[0] += norm2(f[shell]);
788 gmx_sumd(4, buf, cr);
789 ntot = (int)(buf[1] + 0.5);
795 return (ntot ? sqrt(buf[0]/ntot) : 0);
798 static void check_pbc(FILE *fp, rvec x[], int shell)
803 for (m = 0; (m < DIM); m++)
805 if (fabs(x[shell][m]-x[now][m]) > 0.3)
807 pr_rvecs(fp, 0, "SHELL-X", x+now, 5);
813 static void dump_shells(FILE *fp, rvec x[], rvec f[], real ftol, int ns, t_shell s[])
820 for (i = 0; (i < ns); i++)
823 ff2 = iprod(f[shell], f[shell]);
826 fprintf(fp, "SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
827 shell, f[shell][XX], f[shell][YY], f[shell][ZZ], sqrt(ff2));
829 check_pbc(fp, x, shell);
833 static void init_adir(FILE *log, gmx_shellfc_t shfc,
834 gmx_constr_t constr, t_idef *idef, t_inputrec *ir,
835 t_commrec *cr, int dd_ac1,
836 gmx_large_int_t step, t_mdatoms *md, int start, int end,
837 rvec *x_old, rvec *x_init, rvec *x,
838 rvec *f, rvec *acc_dir,
839 gmx_bool bMolPBC, matrix box,
840 real *lambda, real *dvdlambda, t_nrnb *nrnb)
847 unsigned short *ptype;
850 if (DOMAINDECOMP(cr))
858 if (n > shfc->adir_nalloc)
860 shfc->adir_nalloc = over_alloc_dd(n);
861 srenew(shfc->adir_xnold, shfc->adir_nalloc);
862 srenew(shfc->adir_xnew, shfc->adir_nalloc);
864 xnold = shfc->adir_xnold;
865 xnew = shfc->adir_xnew;
871 /* Does NOT work with freeze or acceleration groups (yet) */
872 for (n = start; n < end; n++)
874 w_dt = md->invmass[n]*dt;
876 for (d = 0; d < DIM; d++)
878 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
880 xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
881 xnew[n-start][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
885 xnold[n-start][d] = x[n][d];
886 xnew[n-start][d] = x[n][d];
890 constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
891 x, xnold-start, NULL, bMolPBC, box,
892 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
893 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
894 constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
895 x, xnew-start, NULL, bMolPBC, box,
896 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
897 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
899 for (n = start; n < end; n++)
901 for (d = 0; d < DIM; d++)
904 -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
905 - f[n][d]*md->invmass[n];
907 clear_rvec(acc_dir[n]);
910 /* Project the acceleration on the old bond directions */
911 constrain(log, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, md,
912 x_old, xnew-start, acc_dir, bMolPBC, box,
913 lambda[efptBONDED], &(dvdlambda[efptBONDED]),
914 NULL, NULL, nrnb, econqDeriv_FlexCon, FALSE, 0, 0);
917 int relax_shell_flexcon(FILE *fplog, t_commrec *cr, gmx_bool bVerbose,
918 gmx_large_int_t mdstep, t_inputrec *inputrec,
919 gmx_bool bDoNS, int force_flags,
922 gmx_enerdata_t *enerd, t_fcdata *fcd,
923 t_state *state, rvec f[],
926 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
928 gmx_groups_t *groups,
929 struct gmx_shellfc *shfc,
932 double t, rvec mu_tot,
933 gmx_bool *bConverged,
940 rvec *pos[2], *force[2], *acc_dir = NULL, *x_old = NULL;
944 real ftol, xiH, xiS, dum = 0;
946 gmx_bool bCont, bInit;
947 int nat, dd_ac0, dd_ac1 = 0, i;
948 int start = md->start, homenr = md->homenr, end = start+homenr, cg0, cg1;
949 int nflexcon, g, number_steps, d, Min = 0, count = 0;
950 #define Try (1-Min) /* At start Try = 1 */
952 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
953 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
954 ftol = inputrec->em_tol;
955 number_steps = inputrec->niter;
956 nshell = shfc->nshell;
958 nflexcon = shfc->nflexcon;
962 if (DOMAINDECOMP(cr))
964 nat = dd_natoms_vsite(cr->dd);
967 dd_get_constraint_range(cr->dd, &dd_ac0, &dd_ac1);
968 nat = max(nat, dd_ac1);
976 if (nat > shfc->x_nalloc)
978 /* Allocate local arrays */
979 shfc->x_nalloc = over_alloc_dd(nat);
980 for (i = 0; (i < 2); i++)
982 srenew(shfc->x[i], shfc->x_nalloc);
983 srenew(shfc->f[i], shfc->x_nalloc);
986 for (i = 0; (i < 2); i++)
989 force[i] = shfc->f[i];
992 /* With particle decomposition this code only works
993 * when all particles involved with each shell are in the same cg.
996 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr))
998 /* This is the only time where the coordinates are used
999 * before do_force is called, which normally puts all
1000 * charge groups in the box.
1004 pd_cg_range(cr, &cg0, &cg1);
1011 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, state->box,
1012 &(top->cgs), state->x, fr->cg_cm);
1015 mk_mshift(fplog, graph, fr->ePBC, state->box, state->x);
1019 /* After this all coordinate arrays will contain whole molecules */
1022 shift_self(graph, state->box, state->x);
1027 if (nat > shfc->flex_nalloc)
1029 shfc->flex_nalloc = over_alloc_dd(nat);
1030 srenew(shfc->acc_dir, shfc->flex_nalloc);
1031 srenew(shfc->x_old, shfc->flex_nalloc);
1033 acc_dir = shfc->acc_dir;
1034 x_old = shfc->x_old;
1035 for (i = 0; i < homenr; i++)
1037 for (d = 0; d < DIM; d++)
1040 state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
1045 /* Do a prediction of the shell positions */
1046 if (shfc->bPredict && !bCont)
1048 predict_shells(fplog, state->x, state->v, inputrec->delta_t, nshell, shell,
1049 md->massT, NULL, bInit);
1052 /* do_force expected the charge groups to be in the box */
1055 unshift_self(graph, state->box, state->x);
1058 /* Calculate the forces first time around */
1061 pr_rvecs(debug, 0, "x b4 do_force", state->x + start, homenr);
1063 do_force(fplog, cr, inputrec, mdstep, nrnb, wcycle, top, groups,
1064 state->box, state->x, &state->hist,
1065 force[Min], force_vir, md, enerd, fcd,
1066 state->lambda, graph,
1067 fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1068 (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
1073 init_adir(fplog, shfc,
1074 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1075 shfc->x_old-start, state->x, state->x, force[Min],
1076 shfc->acc_dir-start,
1077 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1079 for (i = start; i < end; i++)
1081 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
1085 Epot[Min] = enerd->term[F_EPOT];
1087 df[Min] = rms_force(cr, shfc->f[Min], nshell, shell, nflexcon, &sf_dir, &Epot[Min]);
1091 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1096 pr_rvecs(debug, 0, "force0", force[Min], md->nr);
1099 if (nshell+nflexcon > 0)
1101 /* Copy x to pos[Min] & pos[Try]: during minimization only the
1102 * shell positions are updated, therefore the other particles must
1105 memcpy(pos[Min], state->x, nat*sizeof(state->x[0]));
1106 memcpy(pos[Try], state->x, nat*sizeof(state->x[0]));
1109 if (bVerbose && MASTER(cr))
1111 print_epot(stdout, mdstep, 0, Epot[Min], df[Min], nflexcon, sf_dir);
1116 fprintf(debug, "%17s: %14.10e\n",
1117 interaction_function[F_EKIN].longname, enerd->term[F_EKIN]);
1118 fprintf(debug, "%17s: %14.10e\n",
1119 interaction_function[F_EPOT].longname, enerd->term[F_EPOT]);
1120 fprintf(debug, "%17s: %14.10e\n",
1121 interaction_function[F_ETOT].longname, enerd->term[F_ETOT]);
1122 fprintf(debug, "SHELLSTEP %s\n", gmx_step_str(mdstep, sbuf));
1125 /* First check whether we should do shells, or whether the force is
1126 * low enough even without minimization.
1128 *bConverged = (df[Min] < ftol);
1130 for (count = 1; (!(*bConverged) && (count < number_steps)); count++)
1134 construct_vsites(vsite, pos[Min], inputrec->delta_t, state->v,
1135 idef->iparams, idef->il,
1136 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1141 init_adir(fplog, shfc,
1142 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1143 x_old-start, state->x, pos[Min], force[Min], acc_dir-start,
1144 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1146 directional_sd(pos[Min], pos[Try], acc_dir-start, start, end,
1150 /* New positions, Steepest descent */
1151 shell_pos_sd(pos[Min], pos[Try], force[Min], nshell, shell, count);
1153 /* do_force expected the charge groups to be in the box */
1156 unshift_self(graph, state->box, pos[Try]);
1161 pr_rvecs(debug, 0, "RELAX: pos[Min] ", pos[Min] + start, homenr);
1162 pr_rvecs(debug, 0, "RELAX: pos[Try] ", pos[Try] + start, homenr);
1164 /* Try the new positions */
1165 do_force(fplog, cr, inputrec, 1, nrnb, wcycle,
1166 top, groups, state->box, pos[Try], &state->hist,
1167 force[Try], force_vir,
1168 md, enerd, fcd, state->lambda, graph,
1169 fr, vsite, mu_tot, t, fp_field, NULL, bBornRadii,
1174 pr_rvecs(debug, 0, "RELAX: force[Min]", force[Min] + start, homenr);
1175 pr_rvecs(debug, 0, "RELAX: force[Try]", force[Try] + start, homenr);
1180 init_adir(fplog, shfc,
1181 constr, idef, inputrec, cr, dd_ac1, mdstep, md, start, end,
1182 x_old-start, state->x, pos[Try], force[Try], acc_dir-start,
1183 fr->bMolPBC, state->box, state->lambda, &dum, nrnb);
1185 for (i = start; i < end; i++)
1187 sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
1191 Epot[Try] = enerd->term[F_EPOT];
1193 df[Try] = rms_force(cr, force[Try], nshell, shell, nflexcon, &sf_dir, &Epot[Try]);
1197 fprintf(debug, "df = %g %g\n", df[Min], df[Try]);
1204 pr_rvecs(debug, 0, "F na do_force", force[Try] + start, homenr);
1208 fprintf(debug, "SHELL ITER %d\n", count);
1209 dump_shells(debug, pos[Try], force[Try], ftol, nshell, shell);
1213 if (bVerbose && MASTER(cr))
1215 print_epot(stdout, mdstep, count, Epot[Try], df[Try], nflexcon, sf_dir);
1218 *bConverged = (df[Try] < ftol);
1220 if ((df[Try] < df[Min]))
1224 fprintf(debug, "Swapping Min and Try\n");
1228 /* Correct the velocities for the flexible constraints */
1229 invdt = 1/inputrec->delta_t;
1230 for (i = start; i < end; i++)
1232 for (d = 0; d < DIM; d++)
1234 state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1242 decrease_step_size(nshell, shell);
1245 if (MASTER(cr) && !(*bConverged))
1247 /* Note that the energies and virial are incorrect when not converged */
1251 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1252 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1255 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1256 gmx_step_str(mdstep, sbuf), number_steps, df[Min]);
1259 /* Copy back the coordinates and the forces */
1260 memcpy(state->x, pos[Min], nat*sizeof(state->x[0]));
1261 memcpy(f, force[Min], nat*sizeof(f[0]));