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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gmx_fatal.h"
60 #include "mtop_util.h"
61 #include "chargegroup.h"
66 atom_id shell; /* The shell id */
67 atom_id nucl1,nucl2,nucl3; /* The nuclei connected to the shell */
68 /* gmx_bool bInterCG; */ /* Coupled to nuclei outside cg? */
69 real k; /* force constant */
70 real k_1; /* 1 over force constant */
76 typedef struct gmx_shellfc {
77 int nshell_gl; /* The number of shells in the system */
78 t_shell *shell_gl; /* All the shells (for DD only) */
79 int *shell_index_gl; /* Global shell index (for DD only) */
80 gmx_bool bInterCG; /* Are there inter charge-group shells? */
81 int nshell; /* The number of local shells */
82 t_shell *shell; /* The local shells */
83 int shell_nalloc; /* The allocation size of shell */
84 gmx_bool bPredict; /* Predict shell positions */
85 gmx_bool bRequireInit; /* Require initialization of shell positions */
86 int nflexcon; /* The number of flexible constraints */
87 rvec *x[2]; /* Array for iterative minimization */
88 rvec *f[2]; /* Array for iterative minimization */
89 int x_nalloc; /* The allocation size of x and f */
90 rvec *acc_dir; /* Acceleration direction for flexcon */
91 rvec *x_old; /* Old coordinates for flexcon */
92 int flex_nalloc; /* The allocation size of acc_dir and x_old */
93 rvec *adir_xnold; /* Work space for init_adir */
94 rvec *adir_xnew; /* Work space for init_adir */
95 int adir_nalloc; /* Work space for init_adir */
99 static void pr_shell(FILE *fplog,int ns,t_shell s[])
103 fprintf(fplog,"SHELL DATA\n");
104 fprintf(fplog,"%5s %8s %5s %5s %5s\n",
105 "Shell","Force k","Nucl1","Nucl2","Nucl3");
106 for(i=0; (i<ns); i++) {
107 fprintf(fplog,"%5d %8.3f %5d",s[i].shell,1.0/s[i].k_1,s[i].nucl1);
109 fprintf(fplog," %5d\n",s[i].nucl2);
110 else if (s[i].nnucl == 3)
111 fprintf(fplog," %5d %5d\n",s[i].nucl2,s[i].nucl3);
117 static void predict_shells(FILE *fplog,rvec x[],rvec v[],real dt,
119 real mass[],gmx_mtop_t *mtop,gmx_bool bInit)
122 real dt_1,dt_2,dt_3,fudge,tm,m1,m2,m3;
124 gmx_mtop_atomlookup_t alook=NULL;
128 alook = gmx_mtop_atomlookup_init(mtop);
131 /* We introduce a fudge factor for performance reasons: with this choice
132 * the initial force on the shells is about a factor of two lower than
139 fprintf(fplog,"RELAX: Using prediction for initial shell placement\n");
148 for(i=0; (i<ns); i++) {
152 switch (s[i].nnucl) {
155 for(m=0; (m<DIM); m++)
156 x[s1][m]+=ptr[n1][m]*dt_1;
165 /* Not the correct masses with FE, but it is just a prediction... */
170 for(m=0; (m<DIM); m++)
171 x[s1][m]+=(m1*ptr[n1][m]+m2*ptr[n2][m])*tm;
182 /* Not the correct masses with FE, but it is just a prediction... */
183 gmx_mtop_atomnr_to_atom(alook,n1,&atom);
185 gmx_mtop_atomnr_to_atom(alook,n2,&atom);
187 gmx_mtop_atomnr_to_atom(alook,n3,&atom);
190 tm = dt_1/(m1+m2+m3);
191 for(m=0; (m<DIM); m++)
192 x[s1][m]+=(m1*ptr[n1][m]+m2*ptr[n2][m]+m3*ptr[n3][m])*tm;
195 gmx_fatal(FARGS,"Shell %d has %d nuclei!",i,s[i].nnucl);
200 gmx_mtop_atomlookup_destroy(alook);
204 gmx_shellfc_t init_shell_flexcon(FILE *fplog,
205 gmx_mtop_t *mtop,int nflexcon,
208 struct gmx_shellfc *shfc;
210 int *shell_index=NULL,*at2cg;
212 int n[eptNR],ns,nshell,nsi;
213 int i,j,nmol,type,mb,mt,a_offset,cg,mol,ftype,nra;
215 int aS,aN=0; /* Shell and nucleus */
216 int bondtypes[] = { F_BONDS, F_HARMONIC, F_CUBICBONDS, F_POLARIZATION, F_ANHARM_POL, F_WATER_POL };
217 #define NBT asize(bondtypes)
219 gmx_mtop_atomloop_block_t aloopb;
220 gmx_mtop_atomloop_all_t aloop;
221 gmx_ffparams_t *ffparams;
222 gmx_molblock_t *molb;
226 /* Count number of shells, and find their indices */
227 for(i=0; (i<eptNR); i++) {
231 aloopb = gmx_mtop_atomloop_block_init(mtop);
232 while (gmx_mtop_atomloop_block_next(aloopb,&atom,&nmol)) {
233 n[atom->ptype] += nmol;
237 /* Print the number of each particle type */
238 for(i=0; (i<eptNR); i++) {
240 fprintf(fplog,"There are: %d %ss\n",n[i],ptype_str[i]);
245 nshell = n[eptShell];
247 if (nshell == 0 && nflexcon == 0) {
252 shfc->nflexcon = nflexcon;
258 /* We have shells: fill the shell data structure */
260 /* Global system sized array, this should be avoided */
261 snew(shell_index,mtop->natoms);
263 aloop = gmx_mtop_atomloop_all_init(mtop);
265 while (gmx_mtop_atomloop_all_next(aloop,&i,&atom)) {
266 if (atom->ptype == eptShell) {
267 shell_index[i] = nshell++;
273 /* Initiate the shell structures */
274 for(i=0; (i<nshell); i++) {
275 shell[i].shell = NO_ATID;
277 shell[i].nucl1 = NO_ATID;
278 shell[i].nucl2 = NO_ATID;
279 shell[i].nucl3 = NO_ATID;
280 /* shell[i].bInterCG=FALSE; */
285 ffparams = &mtop->ffparams;
287 /* Now fill the structures */
288 shfc->bInterCG = FALSE;
291 for(mb=0; mb<mtop->nmolblock; mb++) {
292 molb = &mtop->molblock[mb];
293 molt = &mtop->moltype[molb->type];
296 snew(at2cg,molt->atoms.nr);
297 for(cg=0; cg<cgs->nr; cg++) {
298 for(i=cgs->index[cg]; i<cgs->index[cg+1]; i++) {
303 atom = molt->atoms.atom;
304 for(mol=0; mol<molb->nmol; mol++) {
305 for(j=0; (j<NBT); j++) {
306 ia = molt->ilist[bondtypes[j]].iatoms;
307 for(i=0; (i<molt->ilist[bondtypes[j]].nr); ) {
309 ftype = ffparams->functype[type];
310 nra = interaction_function[ftype].nratoms;
312 /* Check whether we have a bond with a shell */
315 switch (bondtypes[j]) {
321 if (atom[ia[1]].ptype == eptShell) {
325 else if (atom[ia[2]].ptype == eptShell) {
331 aN = ia[4]; /* Dummy */
332 aS = ia[5]; /* Shell */
335 gmx_fatal(FARGS,"Death Horror: %s, %d",__FILE__,__LINE__);
341 /* Check whether one of the particles is a shell... */
342 nsi = shell_index[a_offset+aS];
343 if ((nsi < 0) || (nsi >= nshell))
344 gmx_fatal(FARGS,"nsi is %d should be within 0 - %d. aS = %d",
346 if (shell[nsi].shell == NO_ATID) {
347 shell[nsi].shell = a_offset + aS;
350 else if (shell[nsi].shell != a_offset+aS)
351 gmx_fatal(FARGS,"Weird stuff in %s, %d",__FILE__,__LINE__);
353 if (shell[nsi].nucl1 == NO_ATID) {
354 shell[nsi].nucl1 = a_offset + aN;
355 } else if (shell[nsi].nucl2 == NO_ATID) {
356 shell[nsi].nucl2 = a_offset + aN;
357 } else if (shell[nsi].nucl3 == NO_ATID) {
358 shell[nsi].nucl3 = a_offset + aN;
361 pr_shell(fplog,ns,shell);
362 gmx_fatal(FARGS,"Can not handle more than three bonds per shell\n");
364 if (at2cg[aS] != at2cg[aN]) {
365 /* shell[nsi].bInterCG = TRUE; */
366 shfc->bInterCG = TRUE;
369 switch (bondtypes[j]) {
372 shell[nsi].k += ffparams->iparams[type].harmonic.krA;
375 shell[nsi].k += ffparams->iparams[type].cubic.kb;
379 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
380 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);
381 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/
382 ffparams->iparams[type].polarize.alpha;
385 if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10))
386 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);
387 alpha = (ffparams->iparams[type].wpol.al_x+
388 ffparams->iparams[type].wpol.al_y+
389 ffparams->iparams[type].wpol.al_z)/3.0;
390 shell[nsi].k += sqr(qS)*ONE_4PI_EPS0/alpha;
393 gmx_fatal(FARGS,"Death Horror: %s, %d",__FILE__,__LINE__);
401 a_offset += molt->atoms.nr;
403 /* Done with this molecule type */
407 /* Verify whether it's all correct */
409 gmx_fatal(FARGS,"Something weird with shells. They may not be bonded to something");
411 for(i=0; (i<ns); i++)
412 shell[i].k_1 = 1.0/shell[i].k;
415 pr_shell(debug,ns,shell);
418 shfc->nshell_gl = ns;
419 shfc->shell_gl = shell;
420 shfc->shell_index_gl = shell_index;
422 shfc->bPredict = (getenv("GMX_NOPREDICT") == NULL);
423 shfc->bRequireInit = FALSE;
424 if (!shfc->bPredict) {
426 fprintf(fplog,"\nWill never predict shell positions\n");
428 shfc->bRequireInit = (getenv("GMX_REQUIRE_SHELL_INIT") != NULL);
429 if (shfc->bRequireInit && fplog)
430 fprintf(fplog,"\nWill always initiate shell positions\n");
433 if (shfc->bPredict) {
435 predict_shells(fplog,x,NULL,0,shfc->nshell_gl,shfc->shell_gl,
439 if (shfc->bInterCG) {
441 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");
442 shfc->bPredict = FALSE;
449 void make_local_shells(t_commrec *cr,t_mdatoms *md,
450 struct gmx_shellfc *shfc)
453 int a0,a1,*ind,nshell,i;
454 gmx_domdec_t *dd=NULL;
457 if (DOMAINDECOMP(cr)) {
462 pd_at_range(cr,&a0,&a1);
465 /* Single node: we need all shells, just copy the pointer */
466 shfc->nshell = shfc->nshell_gl;
467 shfc->shell = shfc->shell_gl;
472 ind = shfc->shell_index_gl;
476 for(i=a0; i<a1; i++) {
477 if (md->ptype[i] == eptShell) {
478 if (nshell+1 > shfc->shell_nalloc) {
479 shfc->shell_nalloc = over_alloc_dd(nshell+1);
480 srenew(shell,shfc->shell_nalloc);
483 shell[nshell] = shfc->shell_gl[ind[dd->gatindex[i]]];
485 shell[nshell] = shfc->shell_gl[ind[i]];
487 /* With inter-cg shells we can no do shell prediction,
488 * so we do not need the nuclei numbers.
490 if (!shfc->bInterCG) {
491 shell[nshell].nucl1 = i + shell[nshell].nucl1 - shell[nshell].shell;
492 if (shell[nshell].nnucl > 1)
493 shell[nshell].nucl2 = i + shell[nshell].nucl2 - shell[nshell].shell;
494 if (shell[nshell].nnucl > 2)
495 shell[nshell].nucl3 = i + shell[nshell].nucl3 - shell[nshell].shell;
497 shell[nshell].shell = i;
502 shfc->nshell = nshell;
506 static void do_1pos(rvec xnew,rvec xold,rvec f,real step)
524 static void do_1pos3(rvec xnew,rvec xold,rvec f,rvec step)
542 static void directional_sd(FILE *log,rvec xold[],rvec xnew[],rvec acc_dir[],
543 int start,int homenr,real step)
547 for(i=start; i<homenr; i++)
548 do_1pos(xnew[i],xold[i],acc_dir[i],step);
551 static void shell_pos_sd(FILE *log,rvec xcur[],rvec xnew[],rvec f[],
552 int ns,t_shell s[],int count)
554 const real step_scale_min = 0.8,
555 step_scale_increment = 0.2,
556 step_scale_max = 1.2,
557 step_scale_multiple = (step_scale_max - step_scale_min) / step_scale_increment;
561 real step_min,step_max;
566 for(i=0; (i<ns); i++) {
569 for(d=0; d<DIM; d++) {
570 s[i].step[d] = s[i].k_1;
572 step_min = min(step_min,s[i].step[d]);
573 step_max = max(step_max,s[i].step[d]);
577 for(d=0; d<DIM; d++) {
578 dx = xcur[shell][d] - s[i].xold[d];
579 df = f[shell][d] - s[i].fold[d];
580 /* -dx/df gets used to generate an interpolated value, but would
581 * cause a NaN if df were binary-equal to zero. Values close to
582 * zero won't cause problems (because of the min() and max()), so
583 * just testing for binary inequality is OK. */
587 /* Scale the step size by a factor interpolated from
588 * step_scale_min to step_scale_max, as k_est goes from 0 to
589 * step_scale_multiple * s[i].step[d] */
591 step_scale_min * s[i].step[d] +
592 step_scale_increment * min(step_scale_multiple * s[i].step[d], max(k_est, 0));
597 if (gmx_numzero(dx)) /* 0 == dx */
599 /* Likely this will never happen, but if it does just
600 * don't scale the step. */
604 s[i].step[d] *= step_scale_max;
608 step_min = min(step_min,s[i].step[d]);
609 step_max = max(step_max,s[i].step[d]);
613 copy_rvec(xcur[shell],s[i].xold);
614 copy_rvec(f[shell], s[i].fold);
616 do_1pos3(xnew[shell],xcur[shell],f[shell],s[i].step);
619 fprintf(debug,"shell[%d] = %d\n",i,shell);
620 pr_rvec(debug,0,"fshell",f[shell],DIM,TRUE);
621 pr_rvec(debug,0,"xold",xcur[shell],DIM,TRUE);
622 pr_rvec(debug,0,"step",s[i].step,DIM,TRUE);
623 pr_rvec(debug,0,"xnew",xnew[shell],DIM,TRUE);
627 printf("step %.3e %.3e\n",step_min,step_max);
631 static void decrease_step_size(int nshell,t_shell s[])
635 for(i=0; i<nshell; i++)
636 svmul(0.8,s[i].step,s[i].step);
639 static void print_epot(FILE *fp,gmx_large_int_t mdstep,int count,real epot,real df,
640 int ndir,real sf_dir)
644 fprintf(fp,"MDStep=%5s/%2d EPot: %12.8e, rmsF: %6.2e",
645 gmx_step_str(mdstep,buf),count,epot,df);
647 fprintf(fp,", dir. rmsF: %6.2e\n",sqrt(sf_dir/ndir));
653 static real rms_force(t_commrec *cr,rvec f[],int ns,t_shell s[],
654 int ndir,real *sf_dir,real *Epot)
660 for(i=0; i<ns; i++) {
662 buf[0] += norm2(f[shell]);
671 ntot = (int)(buf[1] + 0.5);
677 return (ntot ? sqrt(buf[0]/ntot) : 0);
680 static void check_pbc(FILE *fp,rvec x[],int shell)
685 for(m=0; (m<DIM); m++)
686 if (fabs(x[shell][m]-x[now][m]) > 0.3) {
687 pr_rvecs(fp,0,"SHELL-X",x+now,5);
692 static void dump_shells(FILE *fp,rvec x[],rvec f[],real ftol,int ns,t_shell s[])
699 for(i=0; (i<ns); i++) {
701 ff2 = iprod(f[shell],f[shell]);
703 fprintf(fp,"SHELL %5d, force %10.5f %10.5f %10.5f, |f| %10.5f\n",
704 shell,f[shell][XX],f[shell][YY],f[shell][ZZ],sqrt(ff2));
705 check_pbc(fp,x,shell);
709 static void init_adir(FILE *log,gmx_shellfc_t shfc,
710 gmx_constr_t constr,t_idef *idef,t_inputrec *ir,
711 t_commrec *cr,int dd_ac1,
712 gmx_large_int_t step,t_mdatoms *md,int start,int end,
713 rvec *x_old,rvec *x_init,rvec *x,
714 rvec *f,rvec *acc_dir,
715 gmx_bool bMolPBC,matrix box,
716 real *lambda,real *dvdlambda,t_nrnb *nrnb)
723 unsigned short *ptype;
726 if (DOMAINDECOMP(cr))
730 if (n > shfc->adir_nalloc) {
731 shfc->adir_nalloc = over_alloc_dd(n);
732 srenew(shfc->adir_xnold,shfc->adir_nalloc);
733 srenew(shfc->adir_xnew ,shfc->adir_nalloc);
735 xnold = shfc->adir_xnold;
736 xnew = shfc->adir_xnew;
742 /* Does NOT work with freeze or acceleration groups (yet) */
743 for (n=start; n<end; n++) {
744 w_dt = md->invmass[n]*dt;
746 for (d=0; d<DIM; d++) {
747 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell)) {
748 xnold[n-start][d] = x[n][d] - (x_init[n][d] - x_old[n][d]);
749 xnew[n-start][d] = 2*x[n][d] - x_old[n][d] + f[n][d]*w_dt*dt;
751 xnold[n-start][d] = x[n][d];
752 xnew[n-start][d] = x[n][d];
756 constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
757 x,xnold-start,NULL,bMolPBC,box,
758 lambda[efptBONDED],&(dvdlambda[efptBONDED]),
759 NULL,NULL,nrnb,econqCoord,FALSE,0,0);
760 constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
761 x,xnew-start,NULL,bMolPBC,box,
762 lambda[efptBONDED],&(dvdlambda[efptBONDED]),
763 NULL,NULL,nrnb,econqCoord,FALSE,0,0);
765 for (n=start; n<end; n++) {
768 -(2*x[n][d]-xnold[n-start][d]-xnew[n-start][d])/sqr(dt)
769 - f[n][d]*md->invmass[n];
770 clear_rvec(acc_dir[n]);
773 /* Project the acceleration on the old bond directions */
774 constrain(log,FALSE,FALSE,constr,idef,ir,NULL,cr,step,0,md,
775 x_old,xnew-start,acc_dir,bMolPBC,box,
776 lambda[efptBONDED],&(dvdlambda[efptBONDED]),
777 NULL,NULL,nrnb,econqDeriv_FlexCon,FALSE,0,0);
780 int relax_shell_flexcon(FILE *fplog,t_commrec *cr,gmx_bool bVerbose,
781 gmx_large_int_t mdstep,t_inputrec *inputrec,
782 gmx_bool bDoNS,int force_flags,
787 gmx_enerdata_t *enerd,t_fcdata *fcd,
788 t_state *state,rvec f[],
791 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
793 gmx_groups_t *groups,
794 struct gmx_shellfc *shfc,
797 double t,rvec mu_tot,
798 int natoms,gmx_bool *bConverged,
805 rvec *pos[2],*force[2],*acc_dir=NULL,*x_old=NULL;
809 real ftol,xiH,xiS,dum=0;
811 gmx_bool bCont,bInit;
812 int nat,dd_ac0,dd_ac1=0,i;
813 int start=md->start,homenr=md->homenr,end=start+homenr,cg0,cg1;
814 int nflexcon,g,number_steps,d,Min=0,count=0;
815 #define Try (1-Min) /* At start Try = 1 */
817 bCont = (mdstep == inputrec->init_step) && inputrec->bContinuation;
818 bInit = (mdstep == inputrec->init_step) || shfc->bRequireInit;
819 ftol = inputrec->em_tol;
820 number_steps = inputrec->niter;
821 nshell = shfc->nshell;
823 nflexcon = shfc->nflexcon;
827 if (DOMAINDECOMP(cr)) {
828 nat = dd_natoms_vsite(cr->dd);
830 dd_get_constraint_range(cr->dd,&dd_ac0,&dd_ac1);
831 nat = max(nat,dd_ac1);
837 if (nat > shfc->x_nalloc) {
838 /* Allocate local arrays */
839 shfc->x_nalloc = over_alloc_dd(nat);
840 for(i=0; (i<2); i++) {
841 srenew(shfc->x[i],shfc->x_nalloc);
842 srenew(shfc->f[i],shfc->x_nalloc);
845 for(i=0; (i<2); i++) {
847 force[i] = shfc->f[i];
850 /* With particle decomposition this code only works
851 * when all particles involved with each shell are in the same cg.
854 if (bDoNS && inputrec->ePBC != epbcNONE && !DOMAINDECOMP(cr)) {
855 /* This is the only time where the coordinates are used
856 * before do_force is called, which normally puts all
857 * charge groups in the box.
859 if (PARTDECOMP(cr)) {
860 pd_cg_range(cr,&cg0,&cg1);
865 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,state->box,
866 &(top->cgs),state->x,fr->cg_cm);
868 mk_mshift(fplog,graph,fr->ePBC,state->box,state->x);
871 /* After this all coordinate arrays will contain whole molecules */
873 shift_self(graph,state->box,state->x);
876 if (nat > shfc->flex_nalloc) {
877 shfc->flex_nalloc = over_alloc_dd(nat);
878 srenew(shfc->acc_dir,shfc->flex_nalloc);
879 srenew(shfc->x_old,shfc->flex_nalloc);
881 acc_dir = shfc->acc_dir;
883 for(i=0; i<homenr; i++) {
886 state->x[start+i][d] - state->v[start+i][d]*inputrec->delta_t;
890 /* Do a prediction of the shell positions */
891 if (shfc->bPredict && !bCont) {
892 predict_shells(fplog,state->x,state->v,inputrec->delta_t,nshell,shell,
893 md->massT,NULL,bInit);
896 /* do_force expected the charge groups to be in the box */
898 unshift_self(graph,state->box,state->x);
900 /* Calculate the forces first time around */
902 pr_rvecs(debug,0,"x b4 do_force",state->x + start,homenr);
904 do_force(fplog,cr,inputrec,mdstep,nrnb,wcycle,top,mtop,groups,
905 state->box,state->x,&state->hist,
906 force[Min],force_vir,md,enerd,fcd,
908 fr,vsite,mu_tot,t,fp_field,NULL,bBornRadii,
909 (bDoNS ? GMX_FORCE_NS : 0) | force_flags);
913 init_adir(fplog,shfc,
914 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
915 shfc->x_old-start,state->x,state->x,force[Min],
917 fr->bMolPBC,state->box,state->lambda,&dum,nrnb);
919 for(i=start; i<end; i++)
920 sf_dir += md->massT[i]*norm2(shfc->acc_dir[i-start]);
923 Epot[Min] = enerd->term[F_EPOT];
925 df[Min]=rms_force(cr,shfc->f[Min],nshell,shell,nflexcon,&sf_dir,&Epot[Min]);
928 fprintf(debug,"df = %g %g\n",df[Min],df[Try]);
932 pr_rvecs(debug,0,"force0",force[Min],md->nr);
935 if (nshell+nflexcon > 0) {
936 /* Copy x to pos[Min] & pos[Try]: during minimization only the
937 * shell positions are updated, therefore the other particles must
940 memcpy(pos[Min],state->x,nat*sizeof(state->x[0]));
941 memcpy(pos[Try],state->x,nat*sizeof(state->x[0]));
944 if (bVerbose && MASTER(cr))
945 print_epot(stdout,mdstep,0,Epot[Min],df[Min],nflexcon,sf_dir);
948 fprintf(debug,"%17s: %14.10e\n",
949 interaction_function[F_EKIN].longname,enerd->term[F_EKIN]);
950 fprintf(debug,"%17s: %14.10e\n",
951 interaction_function[F_EPOT].longname,enerd->term[F_EPOT]);
952 fprintf(debug,"%17s: %14.10e\n",
953 interaction_function[F_ETOT].longname,enerd->term[F_ETOT]);
954 fprintf(debug,"SHELLSTEP %s\n",gmx_step_str(mdstep,sbuf));
957 /* First check whether we should do shells, or whether the force is
958 * low enough even without minimization.
960 *bConverged = (df[Min] < ftol);
962 for(count=1; (!(*bConverged) && (count < number_steps)); count++) {
964 construct_vsites(fplog,vsite,pos[Min],nrnb,inputrec->delta_t,state->v,
965 idef->iparams,idef->il,
966 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
969 init_adir(fplog,shfc,
970 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
971 x_old-start,state->x,pos[Min],force[Min],acc_dir-start,
972 fr->bMolPBC,state->box,state->lambda,&dum,nrnb);
974 directional_sd(fplog,pos[Min],pos[Try],acc_dir-start,start,end,
978 /* New positions, Steepest descent */
979 shell_pos_sd(fplog,pos[Min],pos[Try],force[Min],nshell,shell,count);
981 /* do_force expected the charge groups to be in the box */
983 unshift_self(graph,state->box,pos[Try]);
986 pr_rvecs(debug,0,"RELAX: pos[Min] ",pos[Min] + start,homenr);
987 pr_rvecs(debug,0,"RELAX: pos[Try] ",pos[Try] + start,homenr);
989 /* Try the new positions */
990 do_force(fplog,cr,inputrec,1,nrnb,wcycle,
991 top,mtop,groups,state->box,pos[Try],&state->hist,
992 force[Try],force_vir,
993 md,enerd,fcd,state->lambda,graph,
994 fr,vsite,mu_tot,t,fp_field,NULL,bBornRadii,
998 pr_rvecs(debug,0,"RELAX: force[Min]",force[Min] + start,homenr);
999 pr_rvecs(debug,0,"RELAX: force[Try]",force[Try] + start,homenr);
1003 init_adir(fplog,shfc,
1004 constr,idef,inputrec,cr,dd_ac1,mdstep,md,start,end,
1005 x_old-start,state->x,pos[Try],force[Try],acc_dir-start,
1006 fr->bMolPBC,state->box,state->lambda,&dum,nrnb);
1008 for(i=start; i<end; i++)
1009 sf_dir += md->massT[i]*norm2(acc_dir[i-start]);
1012 Epot[Try] = enerd->term[F_EPOT];
1014 df[Try]=rms_force(cr,force[Try],nshell,shell,nflexcon,&sf_dir,&Epot[Try]);
1017 fprintf(debug,"df = %g %g\n",df[Min],df[Try]);
1021 pr_rvecs(debug,0,"F na do_force",force[Try] + start,homenr);
1023 fprintf(debug,"SHELL ITER %d\n",count);
1024 dump_shells(debug,pos[Try],force[Try],ftol,nshell,shell);
1028 if (bVerbose && MASTER(cr))
1029 print_epot(stdout,mdstep,count,Epot[Try],df[Try],nflexcon,sf_dir);
1031 *bConverged = (df[Try] < ftol);
1033 if ((df[Try] < df[Min])) {
1035 fprintf(debug,"Swapping Min and Try\n");
1037 /* Correct the velocities for the flexible constraints */
1038 invdt = 1/inputrec->delta_t;
1039 for(i=start; i<end; i++) {
1040 for(d=0; d<DIM; d++)
1041 state->v[i][d] += (pos[Try][i][d] - pos[Min][i][d])*invdt;
1046 decrease_step_size(nshell,shell);
1049 if (MASTER(cr) && !(*bConverged)) {
1050 /* Note that the energies and virial are incorrect when not converged */
1053 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1054 gmx_step_str(mdstep,sbuf),number_steps,df[Min]);
1056 "step %s: EM did not converge in %d iterations, RMS force %.3f\n",
1057 gmx_step_str(mdstep,sbuf),number_steps,df[Min]);
1060 /* Copy back the coordinates and the forces */
1061 memcpy(state->x,pos[Min],nat*sizeof(state->x[0]));
1062 memcpy(f,force[Min],nat*sizeof(f[0]));