1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
55 #include "gmx_fatal.h"
72 #include "gmx_wallcycle.h"
73 #include "mtop_util.h"
77 #include "gromacs/linearalgebra/mtxio.h"
78 #include "gromacs/linearalgebra/sparsematrix.h"
89 static em_state_t *init_em_state()
95 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
96 snew(ems->s.lambda,efptNR);
101 static void print_em_start(FILE *fplog,t_commrec *cr,gmx_runtime_t *runtime,
102 gmx_wallcycle_t wcycle,
107 runtime_start(runtime);
109 sprintf(buf,"Started %s",name);
110 print_date_and_time(fplog,cr->nodeid,buf,NULL);
112 wallcycle_start(wcycle,ewcRUN);
114 static void em_time_end(FILE *fplog,t_commrec *cr,gmx_runtime_t *runtime,
115 gmx_wallcycle_t wcycle)
117 wallcycle_stop(wcycle,ewcRUN);
119 runtime_end(runtime);
122 static void sp_header(FILE *out,const char *minimizer,real ftol,int nsteps)
125 fprintf(out,"%s:\n",minimizer);
126 fprintf(out," Tolerance (Fmax) = %12.5e\n",ftol);
127 fprintf(out," Number of steps = %12d\n",nsteps);
130 static void warn_step(FILE *fp,real ftol,gmx_bool bLastStep,gmx_bool bConstrain)
134 fprintf(fp,"\nReached the maximum number of steps before reaching Fmax < %g\n",ftol);
138 fprintf(fp,"\nStepsize too small, or no change in energy.\n"
139 "Converged to machine precision,\n"
140 "but not to the requested precision Fmax < %g\n",
142 if (sizeof(real)<sizeof(double))
144 fprintf(fp,"\nDouble precision normally gives you higher accuracy.\n");
148 fprintf(fp,"You might need to increase your constraint accuracy, or turn\n"
149 "off constraints alltogether (set constraints = none in mdp file)\n");
156 static void print_converged(FILE *fp,const char *alg,real ftol,
157 gmx_large_int_t count,gmx_bool bDone,gmx_large_int_t nsteps,
158 real epot,real fmax, int nfmax, real fnorm)
160 char buf[STEPSTRSIZE];
163 fprintf(fp,"\n%s converged to Fmax < %g in %s steps\n",
164 alg,ftol,gmx_step_str(count,buf));
165 else if(count<nsteps)
166 fprintf(fp,"\n%s converged to machine precision in %s steps,\n"
167 "but did not reach the requested Fmax < %g.\n",
168 alg,gmx_step_str(count,buf),ftol);
170 fprintf(fp,"\n%s did not converge to Fmax < %g in %s steps.\n",
171 alg,ftol,gmx_step_str(count,buf));
174 fprintf(fp,"Potential Energy = %21.14e\n",epot);
175 fprintf(fp,"Maximum force = %21.14e on atom %d\n",fmax,nfmax+1);
176 fprintf(fp,"Norm of force = %21.14e\n",fnorm);
178 fprintf(fp,"Potential Energy = %14.7e\n",epot);
179 fprintf(fp,"Maximum force = %14.7e on atom %d\n",fmax,nfmax+1);
180 fprintf(fp,"Norm of force = %14.7e\n",fnorm);
184 static void get_f_norm_max(t_commrec *cr,
185 t_grpopts *opts,t_mdatoms *mdatoms,rvec *f,
186 real *fnorm,real *fmax,int *a_fmax)
189 real fmax2,fmax2_0,fam;
190 int la_max,a_max,start,end,i,m,gf;
192 /* This routine finds the largest force and returns it.
193 * On parallel machines the global max is taken.
199 start = mdatoms->start;
200 end = mdatoms->homenr + start;
201 if (mdatoms->cFREEZE) {
202 for(i=start; i<end; i++) {
203 gf = mdatoms->cFREEZE[i];
206 if (!opts->nFreeze[gf][m])
215 for(i=start; i<end; i++) {
225 if (la_max >= 0 && DOMAINDECOMP(cr)) {
226 a_max = cr->dd->gatindex[la_max];
231 snew(sum,2*cr->nnodes+1);
232 sum[2*cr->nodeid] = fmax2;
233 sum[2*cr->nodeid+1] = a_max;
234 sum[2*cr->nnodes] = fnorm2;
235 gmx_sumd(2*cr->nnodes+1,sum,cr);
236 fnorm2 = sum[2*cr->nnodes];
237 /* Determine the global maximum */
238 for(i=0; i<cr->nnodes; i++) {
239 if (sum[2*i] > fmax2) {
241 a_max = (int)(sum[2*i+1] + 0.5);
248 *fnorm = sqrt(fnorm2);
255 static void get_state_f_norm_max(t_commrec *cr,
256 t_grpopts *opts,t_mdatoms *mdatoms,
259 get_f_norm_max(cr,opts,mdatoms,ems->f,&ems->fnorm,&ems->fmax,&ems->a_fmax);
262 void init_em(FILE *fplog,const char *title,
263 t_commrec *cr,t_inputrec *ir,
264 t_state *state_global,gmx_mtop_t *top_global,
265 em_state_t *ems,gmx_localtop_t **top,
266 rvec **f,rvec **f_global,
267 t_nrnb *nrnb,rvec mu_tot,
268 t_forcerec *fr,gmx_enerdata_t **enerd,
269 t_graph **graph,t_mdatoms *mdatoms,gmx_global_stat_t *gstat,
270 gmx_vsite_t *vsite,gmx_constr_t constr,
271 int nfile,const t_filenm fnm[],
272 gmx_mdoutf_t **outf,t_mdebin **mdebin)
279 fprintf(fplog,"Initiating %s\n",title);
282 state_global->ngtc = 0;
284 /* Initialize lambda variables */
285 initialize_lambdas(fplog,ir,&(state_global->fep_state),state_global->lambda,NULL);
289 if (DOMAINDECOMP(cr))
291 *top = dd_init_local_top(top_global);
293 dd_init_local_state(cr->dd,state_global,&ems->s);
297 /* Distribute the charge groups over the nodes from the master node */
298 dd_partition_system(fplog,ir->init_step,cr,TRUE,1,
299 state_global,top_global,ir,
300 &ems->s,&ems->f,mdatoms,*top,
301 fr,vsite,NULL,constr,
303 dd_store_state(cr->dd,&ems->s);
307 snew(*f_global,top_global->natoms);
317 snew(*f,top_global->natoms);
319 /* Just copy the state */
320 ems->s = *state_global;
321 snew(ems->s.x,ems->s.nalloc);
322 snew(ems->f,ems->s.nalloc);
323 for(i=0; i<state_global->natoms; i++)
325 copy_rvec(state_global->x[i],ems->s.x[i]);
327 copy_mat(state_global->box,ems->s.box);
329 if (PAR(cr) && ir->eI != eiNM)
331 /* Initialize the particle decomposition and split the topology */
332 *top = split_system(fplog,top_global,ir,cr);
334 pd_cg_range(cr,&fr->cg0,&fr->hcg);
338 *top = gmx_mtop_generate_local_top(top_global,ir);
342 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols)
344 *graph = mk_graph(fplog,&((*top)->idef),0,top_global->natoms,FALSE,FALSE);
353 pd_at_range(cr,&start,&homenr);
359 homenr = top_global->natoms;
361 atoms2md(top_global,ir,0,NULL,start,homenr,mdatoms);
362 update_mdatoms(mdatoms,state_global->lambda[efptFEP]);
366 set_vsite_top(vsite,*top,mdatoms,cr);
372 if (ir->eConstrAlg == econtSHAKE &&
373 gmx_mtop_ftype_count(top_global,F_CONSTR) > 0)
375 gmx_fatal(FARGS,"Can not do energy minimization with %s, use %s\n",
376 econstr_names[econtSHAKE],econstr_names[econtLINCS]);
379 if (!DOMAINDECOMP(cr))
381 set_constraints(constr,*top,ir,mdatoms,cr);
384 if (!ir->bContinuation)
386 /* Constrain the starting coordinates */
388 constrain(PAR(cr) ? NULL : fplog,TRUE,TRUE,constr,&(*top)->idef,
389 ir,NULL,cr,-1,0,mdatoms,
390 ems->s.x,ems->s.x,NULL,ems->s.box,
391 ems->s.lambda[efptFEP],&dvdlambda,
392 NULL,NULL,nrnb,econqCoord,FALSE,0,0);
398 *gstat = global_stat_init(ir);
401 *outf = init_mdoutf(nfile,fnm,0,cr,ir,NULL);
404 init_enerdata(top_global->groups.grps[egcENER].nr,ir->fepvals->n_lambda,
409 /* Init bin for energy stuff */
410 *mdebin = init_mdebin((*outf)->fp_ene,top_global,ir,NULL);
414 calc_shifts(ems->s.box,fr->shift_vec);
417 static void finish_em(FILE *fplog,t_commrec *cr,gmx_mdoutf_t *outf,
418 gmx_runtime_t *runtime,gmx_wallcycle_t wcycle)
420 if (!(cr->duty & DUTY_PME)) {
421 /* Tell the PME only node to finish */
427 em_time_end(fplog,cr,runtime,wcycle);
430 static void swap_em_state(em_state_t *ems1,em_state_t *ems2)
439 static void copy_em_coords(em_state_t *ems,t_state *state)
443 for(i=0; (i<state->natoms); i++)
445 copy_rvec(ems->s.x[i],state->x[i]);
449 static void write_em_traj(FILE *fplog,t_commrec *cr,
451 gmx_bool bX,gmx_bool bF,const char *confout,
452 gmx_mtop_t *top_global,
453 t_inputrec *ir,gmx_large_int_t step,
455 t_state *state_global,rvec *f_global)
459 if ((bX || bF || confout != NULL) && !DOMAINDECOMP(cr))
461 copy_em_coords(state,state_global);
466 if (bX) { mdof_flags |= MDOF_X; }
467 if (bF) { mdof_flags |= MDOF_F; }
468 write_traj(fplog,cr,outf,mdof_flags,
469 top_global,step,(double)step,
470 &state->s,state_global,state->f,f_global,NULL,NULL);
472 if (confout != NULL && MASTER(cr))
474 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
476 /* Make molecules whole only for confout writing */
477 do_pbc_mtop(fplog,ir->ePBC,state_global->box,top_global,
481 write_sto_conf_mtop(confout,
482 *top_global->name,top_global,
483 state_global->x,NULL,ir->ePBC,state_global->box);
487 static void do_em_step(t_commrec *cr,t_inputrec *ir,t_mdatoms *md,
488 em_state_t *ems1,real a,rvec *f,em_state_t *ems2,
489 gmx_constr_t constr,gmx_localtop_t *top,
490 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
491 gmx_large_int_t count)
495 int start,end,gf,i,m;
502 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
503 gmx_incons("state mismatch in do_em_step");
505 s2->flags = s1->flags;
507 if (s2->nalloc != s1->nalloc) {
508 s2->nalloc = s1->nalloc;
509 srenew(s2->x,s1->nalloc);
510 srenew(ems2->f, s1->nalloc);
511 if (s2->flags & (1<<estCGP))
512 srenew(s2->cg_p, s1->nalloc);
515 s2->natoms = s1->natoms;
516 /* Copy free energy state -> is this necessary? */
517 for (i=0;i<efptNR;i++)
519 s2->lambda[i] = s1->lambda[i];
521 copy_mat(s1->box,s2->box);
524 end = md->start + md->homenr;
529 for(i=start; i<end; i++) {
532 for(m=0; m<DIM; m++) {
533 if (ir->opts.nFreeze[gf][m])
536 x2[i][m] = x1[i][m] + a*f[i][m];
540 if (s2->flags & (1<<estCGP)) {
541 /* Copy the CG p vector */
544 for(i=start; i<end; i++)
545 copy_rvec(x1[i],x2[i]);
548 if (DOMAINDECOMP(cr)) {
549 s2->ddp_count = s1->ddp_count;
550 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc) {
551 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
552 srenew(s2->cg_gl,s2->cg_gl_nalloc);
554 s2->ncg_gl = s1->ncg_gl;
555 for(i=0; i<s2->ncg_gl; i++)
556 s2->cg_gl[i] = s1->cg_gl[i];
557 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
561 wallcycle_start(wcycle,ewcCONSTR);
563 constrain(NULL,TRUE,TRUE,constr,&top->idef,
564 ir,NULL,cr,count,0,md,
565 s1->x,s2->x,NULL,s2->box,s2->lambda[efptBONDED],
566 &dvdlambda,NULL,NULL,nrnb,econqCoord,FALSE,0,0);
567 wallcycle_stop(wcycle,ewcCONSTR);
571 static void em_dd_partition_system(FILE *fplog,int step,t_commrec *cr,
572 gmx_mtop_t *top_global,t_inputrec *ir,
573 em_state_t *ems,gmx_localtop_t *top,
574 t_mdatoms *mdatoms,t_forcerec *fr,
575 gmx_vsite_t *vsite,gmx_constr_t constr,
576 t_nrnb *nrnb,gmx_wallcycle_t wcycle)
578 /* Repartition the domain decomposition */
579 wallcycle_start(wcycle,ewcDOMDEC);
580 dd_partition_system(fplog,step,cr,FALSE,1,
583 mdatoms,top,fr,vsite,NULL,constr,
585 dd_store_state(cr->dd,&ems->s);
586 wallcycle_stop(wcycle,ewcDOMDEC);
589 static void evaluate_energy(FILE *fplog,gmx_bool bVerbose,t_commrec *cr,
590 t_state *state_global,gmx_mtop_t *top_global,
591 em_state_t *ems,gmx_localtop_t *top,
592 t_inputrec *inputrec,
593 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
594 gmx_global_stat_t gstat,
595 gmx_vsite_t *vsite,gmx_constr_t constr,
597 t_graph *graph,t_mdatoms *mdatoms,
598 t_forcerec *fr,rvec mu_tot,
599 gmx_enerdata_t *enerd,tensor vir,tensor pres,
600 gmx_large_int_t count,gmx_bool bFirst)
605 tensor force_vir,shake_vir,ekin;
606 real dvdlambda,prescorr,enercorr,dvdlcorr;
609 /* Set the time to the initial time, the time does not change during EM */
610 t = inputrec->init_t;
613 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count)) {
614 /* This the first state or an old state used before the last ns */
618 if (inputrec->nstlist > 0) {
620 } else if (inputrec->nstlist == -1) {
621 nabnsb = natoms_beyond_ns_buffer(inputrec,fr,&top->cgs,NULL,ems->s.x);
623 gmx_sumi(1,&nabnsb,cr);
629 construct_vsites(fplog,vsite,ems->s.x,nrnb,1,NULL,
630 top->idef.iparams,top->idef.il,
631 fr->ePBC,fr->bMolPBC,graph,cr,ems->s.box);
633 if (DOMAINDECOMP(cr)) {
635 /* Repartition the domain decomposition */
636 em_dd_partition_system(fplog,count,cr,top_global,inputrec,
637 ems,top,mdatoms,fr,vsite,constr,
642 /* Calc force & energy on new trial position */
643 /* do_force always puts the charge groups in the box and shifts again
644 * We do not unshift, so molecules are always whole in congrad.c
646 do_force(fplog,cr,inputrec,
647 count,nrnb,wcycle,top,top_global,&top_global->groups,
648 ems->s.box,ems->s.x,&ems->s.hist,
649 ems->f,force_vir,mdatoms,enerd,fcd,
650 ems->s.lambda,graph,fr,vsite,mu_tot,t,NULL,NULL,TRUE,
651 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | GMX_FORCE_VIRIAL |
652 (bNS ? GMX_FORCE_NS | GMX_FORCE_DOLR : 0));
654 /* Clear the unused shake virial and pressure */
655 clear_mat(shake_vir);
658 /* Communicate stuff when parallel */
659 if (PAR(cr) && inputrec->eI != eiNM)
661 wallcycle_start(wcycle,ewcMoveE);
663 global_stat(fplog,gstat,cr,enerd,force_vir,shake_vir,mu_tot,
664 inputrec,NULL,NULL,NULL,1,&terminate,
665 top_global,&ems->s,FALSE,
671 wallcycle_stop(wcycle,ewcMoveE);
674 /* Calculate long range corrections to pressure and energy */
675 calc_dispcorr(fplog,inputrec,fr,count,top_global->natoms,ems->s.box,ems->s.lambda[efptVDW],
676 pres,force_vir,&prescorr,&enercorr,&dvdlcorr);
677 enerd->term[F_DISPCORR] = enercorr;
678 enerd->term[F_EPOT] += enercorr;
679 enerd->term[F_PRES] += prescorr;
680 enerd->term[F_DVDL] += dvdlcorr;
682 ems->epot = enerd->term[F_EPOT];
685 /* Project out the constraint components of the force */
686 wallcycle_start(wcycle,ewcCONSTR);
688 constrain(NULL,FALSE,FALSE,constr,&top->idef,
689 inputrec,NULL,cr,count,0,mdatoms,
690 ems->s.x,ems->f,ems->f,ems->s.box,ems->s.lambda[efptBONDED],&dvdlambda,
691 NULL,&shake_vir,nrnb,econqForceDispl,FALSE,0,0);
692 if (fr->bSepDVDL && fplog)
693 fprintf(fplog,sepdvdlformat,"Constraints",t,dvdlambda);
694 enerd->term[F_DVDL_BONDED] += dvdlambda;
695 m_add(force_vir,shake_vir,vir);
696 wallcycle_stop(wcycle,ewcCONSTR);
698 copy_mat(force_vir,vir);
702 enerd->term[F_PRES] =
703 calc_pres(fr->ePBC,inputrec->nwall,ems->s.box,ekin,vir,pres);
705 sum_dhdl(enerd,ems->s.lambda,inputrec->fepvals);
707 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
709 get_state_f_norm_max(cr,&(inputrec->opts),mdatoms,ems);
713 static double reorder_partsum(t_commrec *cr,t_grpopts *opts,t_mdatoms *mdatoms,
715 em_state_t *s_min,em_state_t *s_b)
719 int ncg,*cg_gl,*index,c,cg,i,a0,a1,a,gf,m;
721 unsigned char *grpnrFREEZE;
724 fprintf(debug,"Doing reorder_partsum\n");
729 cgs_gl = dd_charge_groups_global(cr->dd);
730 index = cgs_gl->index;
732 /* Collect fm in a global vector fmg.
733 * This conflicts with the spirit of domain decomposition,
734 * but to fully optimize this a much more complicated algorithm is required.
736 snew(fmg,mtop->natoms);
738 ncg = s_min->s.ncg_gl;
739 cg_gl = s_min->s.cg_gl;
741 for(c=0; c<ncg; c++) {
745 for(a=a0; a<a1; a++) {
746 copy_rvec(fm[i],fmg[a]);
750 gmx_sum(mtop->natoms*3,fmg[0],cr);
752 /* Now we will determine the part of the sum for the cgs in state s_b */
754 cg_gl = s_b->s.cg_gl;
758 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
759 for(c=0; c<ncg; c++) {
763 for(a=a0; a<a1; a++) {
764 if (mdatoms->cFREEZE && grpnrFREEZE) {
767 for(m=0; m<DIM; m++) {
768 if (!opts->nFreeze[gf][m]) {
769 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
781 static real pr_beta(t_commrec *cr,t_grpopts *opts,t_mdatoms *mdatoms,
783 em_state_t *s_min,em_state_t *s_b)
789 /* This is just the classical Polak-Ribiere calculation of beta;
790 * it looks a bit complicated since we take freeze groups into account,
791 * and might have to sum it in parallel runs.
794 if (!DOMAINDECOMP(cr) ||
795 (s_min->s.ddp_count == cr->dd->ddp_count &&
796 s_b->s.ddp_count == cr->dd->ddp_count)) {
801 /* This part of code can be incorrect with DD,
802 * since the atom ordering in s_b and s_min might differ.
804 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++) {
805 if (mdatoms->cFREEZE)
806 gf = mdatoms->cFREEZE[i];
808 if (!opts->nFreeze[gf][m]) {
809 sum += (fb[i][m] - fm[i][m])*fb[i][m];
813 /* We need to reorder cgs while summing */
814 sum = reorder_partsum(cr,opts,mdatoms,mtop,s_min,s_b);
819 return sum/sqr(s_min->fnorm);
822 double do_cg(FILE *fplog,t_commrec *cr,
823 int nfile,const t_filenm fnm[],
824 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
826 gmx_vsite_t *vsite,gmx_constr_t constr,
828 t_inputrec *inputrec,
829 gmx_mtop_t *top_global,t_fcdata *fcd,
830 t_state *state_global,
832 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
835 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
837 real cpt_period,real max_hours,
838 const char *deviceOptions,
840 gmx_runtime_t *runtime)
842 const char *CG="Polak-Ribiere Conjugate Gradients";
844 em_state_t *s_min,*s_a,*s_b,*s_c;
846 gmx_enerdata_t *enerd;
848 gmx_global_stat_t gstat;
850 rvec *f_global,*p,*sf,*sfm;
851 double gpa,gpb,gpc,tmp,sum[2],minstep;
858 gmx_bool converged,foundlower;
860 gmx_bool do_log=FALSE,do_ene=FALSE,do_x,do_f;
862 int number_steps,neval=0,nstcg=inputrec->nstcgsteep;
864 int i,m,gf,step,nminstep;
869 s_min = init_em_state();
870 s_a = init_em_state();
871 s_b = init_em_state();
872 s_c = init_em_state();
874 /* Init em and store the local state in s_min */
875 init_em(fplog,CG,cr,inputrec,
876 state_global,top_global,s_min,&top,&f,&f_global,
877 nrnb,mu_tot,fr,&enerd,&graph,mdatoms,&gstat,vsite,constr,
878 nfile,fnm,&outf,&mdebin);
880 /* Print to log file */
881 print_em_start(fplog,cr,runtime,wcycle,CG);
883 /* Max number of steps */
884 number_steps=inputrec->nsteps;
887 sp_header(stderr,CG,inputrec->em_tol,number_steps);
889 sp_header(fplog,CG,inputrec->em_tol,number_steps);
891 /* Call the force routine and some auxiliary (neighboursearching etc.) */
892 /* do_force always puts the charge groups in the box and shifts again
893 * We do not unshift, so molecules are always whole in congrad.c
895 evaluate_energy(fplog,bVerbose,cr,
896 state_global,top_global,s_min,top,
897 inputrec,nrnb,wcycle,gstat,
898 vsite,constr,fcd,graph,mdatoms,fr,
899 mu_tot,enerd,vir,pres,-1,TRUE);
903 /* Copy stuff to the energy bin for easy printing etc. */
904 upd_mdebin(mdebin,FALSE,FALSE,(double)step,
905 mdatoms->tmass,enerd,&s_min->s,inputrec->fepvals,inputrec->expandedvals,s_min->s.box,
906 NULL,NULL,vir,pres,NULL,mu_tot,constr);
908 print_ebin_header(fplog,step,step,s_min->s.lambda[efptFEP]);
909 print_ebin(outf->fp_ene,TRUE,FALSE,FALSE,fplog,step,step,eprNORMAL,
910 TRUE,mdebin,fcd,&(top_global->groups),&(inputrec->opts));
914 /* Estimate/guess the initial stepsize */
915 stepsize = inputrec->em_stepsize/s_min->fnorm;
918 fprintf(stderr," F-max = %12.5e on atom %d\n",
919 s_min->fmax,s_min->a_fmax+1);
920 fprintf(stderr," F-Norm = %12.5e\n",
921 s_min->fnorm/sqrt(state_global->natoms));
922 fprintf(stderr,"\n");
923 /* and copy to the log file too... */
924 fprintf(fplog," F-max = %12.5e on atom %d\n",
925 s_min->fmax,s_min->a_fmax+1);
926 fprintf(fplog," F-Norm = %12.5e\n",
927 s_min->fnorm/sqrt(state_global->natoms));
930 /* Start the loop over CG steps.
931 * Each successful step is counted, and we continue until
932 * we either converge or reach the max number of steps.
935 for(step=0; (number_steps<0 || (number_steps>=0 && step<=number_steps)) && !converged;step++) {
937 /* start taking steps in a new direction
938 * First time we enter the routine, beta=0, and the direction is
939 * simply the negative gradient.
942 /* Calculate the new direction in p, and the gradient in this direction, gpa */
947 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++) {
948 if (mdatoms->cFREEZE)
949 gf = mdatoms->cFREEZE[i];
950 for(m=0; m<DIM; m++) {
951 if (!inputrec->opts.nFreeze[gf][m]) {
952 p[i][m] = sf[i][m] + beta*p[i][m];
953 gpa -= p[i][m]*sf[i][m];
954 /* f is negative gradient, thus the sign */
961 /* Sum the gradient along the line across CPUs */
965 /* Calculate the norm of the search vector */
966 get_f_norm_max(cr,&(inputrec->opts),mdatoms,p,&pnorm,NULL,NULL);
968 /* Just in case stepsize reaches zero due to numerical precision... */
970 stepsize = inputrec->em_stepsize/pnorm;
973 * Double check the value of the derivative in the search direction.
974 * If it is positive it must be due to the old information in the
975 * CG formula, so just remove that and start over with beta=0.
976 * This corresponds to a steepest descent step.
980 step--; /* Don't count this step since we are restarting */
981 continue; /* Go back to the beginning of the big for-loop */
984 /* Calculate minimum allowed stepsize, before the average (norm)
985 * relative change in coordinate is smaller than precision
988 for (i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++) {
989 for(m=0; m<DIM; m++) {
990 tmp = fabs(s_min->s.x[i][m]);
997 /* Add up from all CPUs */
999 gmx_sumd(1,&minstep,cr);
1001 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1003 if(stepsize<minstep) {
1008 /* Write coordinates if necessary */
1009 do_x = do_per_step(step,inputrec->nstxout);
1010 do_f = do_per_step(step,inputrec->nstfout);
1012 write_em_traj(fplog,cr,outf,do_x,do_f,NULL,
1013 top_global,inputrec,step,
1014 s_min,state_global,f_global);
1016 /* Take a step downhill.
1017 * In theory, we should minimize the function along this direction.
1018 * That is quite possible, but it turns out to take 5-10 function evaluations
1019 * for each line. However, we dont really need to find the exact minimum -
1020 * it is much better to start a new CG step in a modified direction as soon
1021 * as we are close to it. This will save a lot of energy evaluations.
1023 * In practice, we just try to take a single step.
1024 * If it worked (i.e. lowered the energy), we increase the stepsize but
1025 * the continue straight to the next CG step without trying to find any minimum.
1026 * If it didn't work (higher energy), there must be a minimum somewhere between
1027 * the old position and the new one.
1029 * Due to the finite numerical accuracy, it turns out that it is a good idea
1030 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1031 * This leads to lower final energies in the tests I've done. / Erik
1033 s_a->epot = s_min->epot;
1035 c = a + stepsize; /* reference position along line is zero */
1037 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count) {
1038 em_dd_partition_system(fplog,step,cr,top_global,inputrec,
1039 s_min,top,mdatoms,fr,vsite,constr,
1043 /* Take a trial step (new coords in s_c) */
1044 do_em_step(cr,inputrec,mdatoms,s_min,c,s_min->s.cg_p,s_c,
1045 constr,top,nrnb,wcycle,-1);
1048 /* Calculate energy for the trial step */
1049 evaluate_energy(fplog,bVerbose,cr,
1050 state_global,top_global,s_c,top,
1051 inputrec,nrnb,wcycle,gstat,
1052 vsite,constr,fcd,graph,mdatoms,fr,
1053 mu_tot,enerd,vir,pres,-1,FALSE);
1055 /* Calc derivative along line */
1059 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++) {
1060 for(m=0; m<DIM; m++)
1061 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1063 /* Sum the gradient along the line across CPUs */
1065 gmx_sumd(1,&gpc,cr);
1067 /* This is the max amount of increase in energy we tolerate */
1068 tmp=sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1070 /* Accept the step if the energy is lower, or if it is not significantly higher
1071 * and the line derivative is still negative.
1073 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp))) {
1075 /* Great, we found a better energy. Increase step for next iteration
1076 * if we are still going down, decrease it otherwise
1079 stepsize *= 1.618034; /* The golden section */
1081 stepsize *= 0.618034; /* 1/golden section */
1083 /* New energy is the same or higher. We will have to do some work
1084 * to find a smaller value in the interval. Take smaller step next time!
1087 stepsize *= 0.618034;
1093 /* OK, if we didn't find a lower value we will have to locate one now - there must
1094 * be one in the interval [a=0,c].
1095 * The same thing is valid here, though: Don't spend dozens of iterations to find
1096 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1097 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1099 * I also have a safeguard for potentially really patological functions so we never
1100 * take more than 20 steps before we give up ...
1102 * If we already found a lower value we just skip this step and continue to the update.
1108 /* Select a new trial point.
1109 * If the derivatives at points a & c have different sign we interpolate to zero,
1110 * otherwise just do a bisection.
1113 b = a + gpa*(a-c)/(gpc-gpa);
1117 /* safeguard if interpolation close to machine accuracy causes errors:
1118 * never go outside the interval
1123 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count) {
1124 /* Reload the old state */
1125 em_dd_partition_system(fplog,-1,cr,top_global,inputrec,
1126 s_min,top,mdatoms,fr,vsite,constr,
1130 /* Take a trial step to this new point - new coords in s_b */
1131 do_em_step(cr,inputrec,mdatoms,s_min,b,s_min->s.cg_p,s_b,
1132 constr,top,nrnb,wcycle,-1);
1135 /* Calculate energy for the trial step */
1136 evaluate_energy(fplog,bVerbose,cr,
1137 state_global,top_global,s_b,top,
1138 inputrec,nrnb,wcycle,gstat,
1139 vsite,constr,fcd,graph,mdatoms,fr,
1140 mu_tot,enerd,vir,pres,-1,FALSE);
1142 /* p does not change within a step, but since the domain decomposition
1143 * might change, we have to use cg_p of s_b here.
1148 for(i=mdatoms->start; i<mdatoms->start+mdatoms->homenr; i++) {
1149 for(m=0; m<DIM; m++)
1150 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1152 /* Sum the gradient along the line across CPUs */
1154 gmx_sumd(1,&gpb,cr);
1157 fprintf(debug,"CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1158 s_a->epot,s_b->epot,s_c->epot,gpb);
1160 epot_repl = s_b->epot;
1162 /* Keep one of the intervals based on the value of the derivative at the new point */
1164 /* Replace c endpoint with b */
1165 swap_em_state(s_b,s_c);
1169 /* Replace a endpoint with b */
1170 swap_em_state(s_b,s_a);
1176 * Stop search as soon as we find a value smaller than the endpoints.
1177 * Never run more than 20 steps, no matter what.
1180 } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1183 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1185 /* OK. We couldn't find a significantly lower energy.
1186 * If beta==0 this was steepest descent, and then we give up.
1187 * If not, set beta=0 and restart with steepest descent before quitting.
1194 /* Reset memory before giving up */
1200 /* Select min energy state of A & C, put the best in B.
1202 if (s_c->epot < s_a->epot) {
1204 fprintf(debug,"CGE: C (%f) is lower than A (%f), moving C to B\n",
1205 s_c->epot,s_a->epot);
1206 swap_em_state(s_b,s_c);
1211 fprintf(debug,"CGE: A (%f) is lower than C (%f), moving A to B\n",
1212 s_a->epot,s_c->epot);
1213 swap_em_state(s_b,s_a);
1220 fprintf(debug,"CGE: Found a lower energy %f, moving C to B\n",
1222 swap_em_state(s_b,s_c);
1227 /* new search direction */
1228 /* beta = 0 means forget all memory and restart with steepest descents. */
1229 if (nstcg && ((step % nstcg)==0))
1232 /* s_min->fnorm cannot be zero, because then we would have converged
1236 /* Polak-Ribiere update.
1237 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1239 beta = pr_beta(cr,&inputrec->opts,mdatoms,top_global,s_min,s_b);
1241 /* Limit beta to prevent oscillations */
1242 if (fabs(beta) > 5.0)
1246 /* update positions */
1247 swap_em_state(s_min,s_b);
1250 /* Print it if necessary */
1253 fprintf(stderr,"\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1254 step,s_min->epot,s_min->fnorm/sqrt(state_global->natoms),
1255 s_min->fmax,s_min->a_fmax+1);
1256 /* Store the new (lower) energies */
1257 upd_mdebin(mdebin,FALSE,FALSE,(double)step,
1258 mdatoms->tmass,enerd,&s_min->s,inputrec->fepvals,inputrec->expandedvals,s_min->s.box,
1259 NULL,NULL,vir,pres,NULL,mu_tot,constr);
1261 do_log = do_per_step(step,inputrec->nstlog);
1262 do_ene = do_per_step(step,inputrec->nstenergy);
1264 print_ebin_header(fplog,step,step,s_min->s.lambda[efptFEP]);
1265 print_ebin(outf->fp_ene,do_ene,FALSE,FALSE,
1266 do_log ? fplog : NULL,step,step,eprNORMAL,
1267 TRUE,mdebin,fcd,&(top_global->groups),&(inputrec->opts));
1270 /* Stop when the maximum force lies below tolerance.
1271 * If we have reached machine precision, converged is already set to true.
1273 converged = converged || (s_min->fmax < inputrec->em_tol);
1275 } /* End of the loop */
1278 step--; /* we never took that last step in this case */
1280 if (s_min->fmax > inputrec->em_tol)
1284 warn_step(stderr,inputrec->em_tol,step-1==number_steps,FALSE);
1285 warn_step(fplog ,inputrec->em_tol,step-1==number_steps,FALSE);
1291 /* If we printed energy and/or logfile last step (which was the last step)
1292 * we don't have to do it again, but otherwise print the final values.
1295 /* Write final value to log since we didn't do anything the last step */
1296 print_ebin_header(fplog,step,step,s_min->s.lambda[efptFEP]);
1298 if (!do_ene || !do_log) {
1299 /* Write final energy file entries */
1300 print_ebin(outf->fp_ene,!do_ene,FALSE,FALSE,
1301 !do_log ? fplog : NULL,step,step,eprNORMAL,
1302 TRUE,mdebin,fcd,&(top_global->groups),&(inputrec->opts));
1306 /* Print some stuff... */
1308 fprintf(stderr,"\nwriting lowest energy coordinates.\n");
1311 * For accurate normal mode calculation it is imperative that we
1312 * store the last conformation into the full precision binary trajectory.
1314 * However, we should only do it if we did NOT already write this step
1315 * above (which we did if do_x or do_f was true).
1317 do_x = !do_per_step(step,inputrec->nstxout);
1318 do_f = (inputrec->nstfout > 0 && !do_per_step(step,inputrec->nstfout));
1320 write_em_traj(fplog,cr,outf,do_x,do_f,ftp2fn(efSTO,nfile,fnm),
1321 top_global,inputrec,step,
1322 s_min,state_global,f_global);
1324 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1327 print_converged(stderr,CG,inputrec->em_tol,step,converged,number_steps,
1328 s_min->epot,s_min->fmax,s_min->a_fmax,fnormn);
1329 print_converged(fplog,CG,inputrec->em_tol,step,converged,number_steps,
1330 s_min->epot,s_min->fmax,s_min->a_fmax,fnormn);
1332 fprintf(fplog,"\nPerformed %d energy evaluations in total.\n",neval);
1335 finish_em(fplog,cr,outf,runtime,wcycle);
1337 /* To print the actual number of steps we needed somewhere */
1338 runtime->nsteps_done = step;
1341 } /* That's all folks */
1344 double do_lbfgs(FILE *fplog,t_commrec *cr,
1345 int nfile,const t_filenm fnm[],
1346 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
1348 gmx_vsite_t *vsite,gmx_constr_t constr,
1350 t_inputrec *inputrec,
1351 gmx_mtop_t *top_global,t_fcdata *fcd,
1354 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1357 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1358 gmx_membed_t membed,
1359 real cpt_period,real max_hours,
1360 const char *deviceOptions,
1361 unsigned long Flags,
1362 gmx_runtime_t *runtime)
1364 static const char *LBFGS="Low-Memory BFGS Minimizer";
1366 gmx_localtop_t *top;
1367 gmx_enerdata_t *enerd;
1369 gmx_global_stat_t gstat;
1372 int ncorr,nmaxcorr,point,cp,neval,nminstep;
1373 double stepsize,gpa,gpb,gpc,tmp,minstep;
1374 real *rho,*alpha,*ff,*xx,*p,*s,*lastx,*lastf,**dx,**dg;
1375 real *xa,*xb,*xc,*fa,*fb,*fc,*xtmp,*ftmp;
1376 real a,b,c,maxdelta,delta;
1377 real diag,Epot0,Epot,EpotA,EpotB,EpotC;
1378 real dgdx,dgdg,sq,yr,beta;
1380 gmx_bool converged,first;
1383 gmx_bool do_log,do_ene,do_x,do_f,foundlower,*frozen;
1385 int start,end,number_steps;
1387 int i,k,m,n,nfmax,gf,step;
1393 gmx_fatal(FARGS,"Cannot do parallel L-BFGS Minimization - yet.\n");
1395 n = 3*state->natoms;
1396 nmaxcorr = inputrec->nbfgscorr;
1398 /* Allocate memory */
1399 /* Use pointers to real so we dont have to loop over both atoms and
1400 * dimensions all the time...
1401 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1402 * that point to the same memory.
1416 snew(alpha,nmaxcorr);
1419 for(i=0;i<nmaxcorr;i++)
1423 for(i=0;i<nmaxcorr;i++)
1430 init_em(fplog,LBFGS,cr,inputrec,
1431 state,top_global,&ems,&top,&f,&f_global,
1432 nrnb,mu_tot,fr,&enerd,&graph,mdatoms,&gstat,vsite,constr,
1433 nfile,fnm,&outf,&mdebin);
1434 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1435 * so we free some memory again.
1440 xx = (real *)state->x;
1443 start = mdatoms->start;
1444 end = mdatoms->homenr + start;
1446 /* Print to log file */
1447 print_em_start(fplog,cr,runtime,wcycle,LBFGS);
1449 do_log = do_ene = do_x = do_f = TRUE;
1451 /* Max number of steps */
1452 number_steps=inputrec->nsteps;
1454 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1456 for(i=start; i<end; i++) {
1457 if (mdatoms->cFREEZE)
1458 gf = mdatoms->cFREEZE[i];
1459 for(m=0; m<DIM; m++)
1460 frozen[3*i+m]=inputrec->opts.nFreeze[gf][m];
1463 sp_header(stderr,LBFGS,inputrec->em_tol,number_steps);
1465 sp_header(fplog,LBFGS,inputrec->em_tol,number_steps);
1468 construct_vsites(fplog,vsite,state->x,nrnb,1,NULL,
1469 top->idef.iparams,top->idef.il,
1470 fr->ePBC,fr->bMolPBC,graph,cr,state->box);
1472 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1473 /* do_force always puts the charge groups in the box and shifts again
1474 * We do not unshift, so molecules are always whole
1479 evaluate_energy(fplog,bVerbose,cr,
1480 state,top_global,&ems,top,
1481 inputrec,nrnb,wcycle,gstat,
1482 vsite,constr,fcd,graph,mdatoms,fr,
1483 mu_tot,enerd,vir,pres,-1,TRUE);
1487 /* Copy stuff to the energy bin for easy printing etc. */
1488 upd_mdebin(mdebin,FALSE,FALSE,(double)step,
1489 mdatoms->tmass,enerd,state,inputrec->fepvals,inputrec->expandedvals,state->box,
1490 NULL,NULL,vir,pres,NULL,mu_tot,constr);
1492 print_ebin_header(fplog,step,step,state->lambda[efptFEP]);
1493 print_ebin(outf->fp_ene,TRUE,FALSE,FALSE,fplog,step,step,eprNORMAL,
1494 TRUE,mdebin,fcd,&(top_global->groups),&(inputrec->opts));
1498 /* This is the starting energy */
1499 Epot = enerd->term[F_EPOT];
1505 /* Set the initial step.
1506 * since it will be multiplied by the non-normalized search direction
1507 * vector (force vector the first time), we scale it by the
1508 * norm of the force.
1512 fprintf(stderr,"Using %d BFGS correction steps.\n\n",nmaxcorr);
1513 fprintf(stderr," F-max = %12.5e on atom %d\n",fmax,nfmax+1);
1514 fprintf(stderr," F-Norm = %12.5e\n",fnorm/sqrt(state->natoms));
1515 fprintf(stderr,"\n");
1516 /* and copy to the log file too... */
1517 fprintf(fplog,"Using %d BFGS correction steps.\n\n",nmaxcorr);
1518 fprintf(fplog," F-max = %12.5e on atom %d\n",fmax,nfmax+1);
1519 fprintf(fplog," F-Norm = %12.5e\n",fnorm/sqrt(state->natoms));
1520 fprintf(fplog,"\n");
1526 dx[point][i] = ff[i]; /* Initial search direction */
1530 stepsize = 1.0/fnorm;
1533 /* Start the loop over BFGS steps.
1534 * Each successful step is counted, and we continue until
1535 * we either converge or reach the max number of steps.
1540 /* Set the gradient from the force */
1542 for(step=0; (number_steps<0 || (number_steps>=0 && step<=number_steps)) && !converged; step++) {
1544 /* Write coordinates if necessary */
1545 do_x = do_per_step(step,inputrec->nstxout);
1546 do_f = do_per_step(step,inputrec->nstfout);
1551 mdof_flags |= MDOF_X;
1556 mdof_flags |= MDOF_F;
1559 write_traj(fplog,cr,outf,mdof_flags,
1560 top_global,step,(real)step,state,state,f,f,NULL,NULL);
1562 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1564 /* pointer to current direction - point=0 first time here */
1567 /* calculate line gradient */
1568 for(gpa=0,i=0;i<n;i++)
1571 /* Calculate minimum allowed stepsize, before the average (norm)
1572 * relative change in coordinate is smaller than precision
1574 for(minstep=0,i=0;i<n;i++) {
1581 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1583 if(stepsize<minstep) {
1588 /* Store old forces and coordinates */
1600 /* Take a step downhill.
1601 * In theory, we should minimize the function along this direction.
1602 * That is quite possible, but it turns out to take 5-10 function evaluations
1603 * for each line. However, we dont really need to find the exact minimum -
1604 * it is much better to start a new BFGS step in a modified direction as soon
1605 * as we are close to it. This will save a lot of energy evaluations.
1607 * In practice, we just try to take a single step.
1608 * If it worked (i.e. lowered the energy), we increase the stepsize but
1609 * the continue straight to the next BFGS step without trying to find any minimum.
1610 * If it didn't work (higher energy), there must be a minimum somewhere between
1611 * the old position and the new one.
1613 * Due to the finite numerical accuracy, it turns out that it is a good idea
1614 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1615 * This leads to lower final energies in the tests I've done. / Erik
1620 c = a + stepsize; /* reference position along line is zero */
1622 /* Check stepsize first. We do not allow displacements
1623 * larger than emstep.
1633 if(maxdelta>inputrec->em_stepsize)
1635 } while(maxdelta>inputrec->em_stepsize);
1637 /* Take a trial step */
1639 xc[i] = lastx[i] + c*s[i];
1642 /* Calculate energy for the trial step */
1643 ems.s.x = (rvec *)xc;
1645 evaluate_energy(fplog,bVerbose,cr,
1646 state,top_global,&ems,top,
1647 inputrec,nrnb,wcycle,gstat,
1648 vsite,constr,fcd,graph,mdatoms,fr,
1649 mu_tot,enerd,vir,pres,step,FALSE);
1652 /* Calc derivative along line */
1653 for(gpc=0,i=0; i<n; i++) {
1654 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1656 /* Sum the gradient along the line across CPUs */
1658 gmx_sumd(1,&gpc,cr);
1660 /* This is the max amount of increase in energy we tolerate */
1661 tmp=sqrt(GMX_REAL_EPS)*fabs(EpotA);
1663 /* Accept the step if the energy is lower, or if it is not significantly higher
1664 * and the line derivative is still negative.
1666 if(EpotC<EpotA || (gpc<0 && EpotC<(EpotA+tmp))) {
1668 /* Great, we found a better energy. Increase step for next iteration
1669 * if we are still going down, decrease it otherwise
1672 stepsize *= 1.618034; /* The golden section */
1674 stepsize *= 0.618034; /* 1/golden section */
1676 /* New energy is the same or higher. We will have to do some work
1677 * to find a smaller value in the interval. Take smaller step next time!
1680 stepsize *= 0.618034;
1683 /* OK, if we didn't find a lower value we will have to locate one now - there must
1684 * be one in the interval [a=0,c].
1685 * The same thing is valid here, though: Don't spend dozens of iterations to find
1686 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1687 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1689 * I also have a safeguard for potentially really patological functions so we never
1690 * take more than 20 steps before we give up ...
1692 * If we already found a lower value we just skip this step and continue to the update.
1699 /* Select a new trial point.
1700 * If the derivatives at points a & c have different sign we interpolate to zero,
1701 * otherwise just do a bisection.
1705 b = a + gpa*(a-c)/(gpc-gpa);
1709 /* safeguard if interpolation close to machine accuracy causes errors:
1710 * never go outside the interval
1715 /* Take a trial step */
1717 xb[i] = lastx[i] + b*s[i];
1720 /* Calculate energy for the trial step */
1721 ems.s.x = (rvec *)xb;
1723 evaluate_energy(fplog,bVerbose,cr,
1724 state,top_global,&ems,top,
1725 inputrec,nrnb,wcycle,gstat,
1726 vsite,constr,fcd,graph,mdatoms,fr,
1727 mu_tot,enerd,vir,pres,step,FALSE);
1732 for(gpb=0,i=0; i<n; i++)
1733 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
1735 /* Sum the gradient along the line across CPUs */
1737 gmx_sumd(1,&gpb,cr);
1739 /* Keep one of the intervals based on the value of the derivative at the new point */
1741 /* Replace c endpoint with b */
1745 /* swap coord pointers b/c */
1753 /* Replace a endpoint with b */
1757 /* swap coord pointers a/b */
1767 * Stop search as soon as we find a value smaller than the endpoints,
1768 * or if the tolerance is below machine precision.
1769 * Never run more than 20 steps, no matter what.
1772 } while((EpotB>EpotA || EpotB>EpotC) && (nminstep<20));
1774 if(fabs(EpotB-Epot0)<GMX_REAL_EPS || nminstep>=20) {
1775 /* OK. We couldn't find a significantly lower energy.
1776 * If ncorr==0 this was steepest descent, and then we give up.
1777 * If not, reset memory to restart as steepest descent before quitting.
1786 /* Search in gradient direction */
1789 /* Reset stepsize */
1790 stepsize = 1.0/fnorm;
1795 /* Select min energy state of A & C, put the best in xx/ff/Epot
1826 /* Update the memory information, and calculate a new
1827 * approximation of the inverse hessian
1830 /* Have new data in Epot, xx, ff */
1835 dg[point][i]=lastf[i]-ff[i];
1836 dx[point][i]*=stepsize;
1842 dgdg+=dg[point][i]*dg[point][i];
1843 dgdx+=dg[point][i]*dx[point][i];
1848 rho[point]=1.0/dgdx;
1860 /* Recursive update. First go back over the memory points */
1861 for(k=0;k<ncorr;k++) {
1870 alpha[cp]=rho[cp]*sq;
1873 p[i] -= alpha[cp]*dg[cp][i];
1879 /* And then go forward again */
1880 for(k=0;k<ncorr;k++) {
1883 yr += p[i]*dg[cp][i];
1886 beta = alpha[cp]-beta;
1889 p[i] += beta*dx[cp][i];
1898 dx[point][i] = p[i];
1904 /* Test whether the convergence criterion is met */
1905 get_f_norm_max(cr,&(inputrec->opts),mdatoms,f,&fnorm,&fmax,&nfmax);
1907 /* Print it if necessary */
1910 fprintf(stderr,"\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1911 step,Epot,fnorm/sqrt(state->natoms),fmax,nfmax+1);
1912 /* Store the new (lower) energies */
1913 upd_mdebin(mdebin,FALSE,FALSE,(double)step,
1914 mdatoms->tmass,enerd,state,inputrec->fepvals,inputrec->expandedvals,state->box,
1915 NULL,NULL,vir,pres,NULL,mu_tot,constr);
1916 do_log = do_per_step(step,inputrec->nstlog);
1917 do_ene = do_per_step(step,inputrec->nstenergy);
1919 print_ebin_header(fplog,step,step,state->lambda[efptFEP]);
1920 print_ebin(outf->fp_ene,do_ene,FALSE,FALSE,
1921 do_log ? fplog : NULL,step,step,eprNORMAL,
1922 TRUE,mdebin,fcd,&(top_global->groups),&(inputrec->opts));
1925 /* Stop when the maximum force lies below tolerance.
1926 * If we have reached machine precision, converged is already set to true.
1929 converged = converged || (fmax < inputrec->em_tol);
1931 } /* End of the loop */
1934 step--; /* we never took that last step in this case */
1936 if(fmax>inputrec->em_tol)
1940 warn_step(stderr,inputrec->em_tol,step-1==number_steps,FALSE);
1941 warn_step(fplog ,inputrec->em_tol,step-1==number_steps,FALSE);
1946 /* If we printed energy and/or logfile last step (which was the last step)
1947 * we don't have to do it again, but otherwise print the final values.
1949 if(!do_log) /* Write final value to log since we didn't do anythin last step */
1950 print_ebin_header(fplog,step,step,state->lambda[efptFEP]);
1951 if(!do_ene || !do_log) /* Write final energy file entries */
1952 print_ebin(outf->fp_ene,!do_ene,FALSE,FALSE,
1953 !do_log ? fplog : NULL,step,step,eprNORMAL,
1954 TRUE,mdebin,fcd,&(top_global->groups),&(inputrec->opts));
1956 /* Print some stuff... */
1958 fprintf(stderr,"\nwriting lowest energy coordinates.\n");
1961 * For accurate normal mode calculation it is imperative that we
1962 * store the last conformation into the full precision binary trajectory.
1964 * However, we should only do it if we did NOT already write this step
1965 * above (which we did if do_x or do_f was true).
1967 do_x = !do_per_step(step,inputrec->nstxout);
1968 do_f = !do_per_step(step,inputrec->nstfout);
1969 write_em_traj(fplog,cr,outf,do_x,do_f,ftp2fn(efSTO,nfile,fnm),
1970 top_global,inputrec,step,
1974 print_converged(stderr,LBFGS,inputrec->em_tol,step,converged,
1975 number_steps,Epot,fmax,nfmax,fnorm/sqrt(state->natoms));
1976 print_converged(fplog,LBFGS,inputrec->em_tol,step,converged,
1977 number_steps,Epot,fmax,nfmax,fnorm/sqrt(state->natoms));
1979 fprintf(fplog,"\nPerformed %d energy evaluations in total.\n",neval);
1982 finish_em(fplog,cr,outf,runtime,wcycle);
1984 /* To print the actual number of steps we needed somewhere */
1985 runtime->nsteps_done = step;
1988 } /* That's all folks */
1991 double do_steep(FILE *fplog,t_commrec *cr,
1992 int nfile, const t_filenm fnm[],
1993 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
1995 gmx_vsite_t *vsite,gmx_constr_t constr,
1997 t_inputrec *inputrec,
1998 gmx_mtop_t *top_global,t_fcdata *fcd,
1999 t_state *state_global,
2001 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
2004 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2005 gmx_membed_t membed,
2006 real cpt_period,real max_hours,
2007 const char *deviceOptions,
2008 unsigned long Flags,
2009 gmx_runtime_t *runtime)
2011 const char *SD="Steepest Descents";
2012 em_state_t *s_min,*s_try;
2014 gmx_localtop_t *top;
2015 gmx_enerdata_t *enerd;
2017 gmx_global_stat_t gstat;
2019 real stepsize,constepsize;
2020 real ustep,dvdlambda,fnormn;
2023 gmx_bool bDone,bAbort,do_x,do_f;
2028 int steps_accepted=0;
2032 s_min = init_em_state();
2033 s_try = init_em_state();
2035 /* Init em and store the local state in s_try */
2036 init_em(fplog,SD,cr,inputrec,
2037 state_global,top_global,s_try,&top,&f,&f_global,
2038 nrnb,mu_tot,fr,&enerd,&graph,mdatoms,&gstat,vsite,constr,
2039 nfile,fnm,&outf,&mdebin);
2041 /* Print to log file */
2042 print_em_start(fplog,cr,runtime,wcycle,SD);
2044 /* Set variables for stepsize (in nm). This is the largest
2045 * step that we are going to make in any direction.
2047 ustep = inputrec->em_stepsize;
2050 /* Max number of steps */
2051 nsteps = inputrec->nsteps;
2054 /* Print to the screen */
2055 sp_header(stderr,SD,inputrec->em_tol,nsteps);
2057 sp_header(fplog,SD,inputrec->em_tol,nsteps);
2059 /**** HERE STARTS THE LOOP ****
2060 * count is the counter for the number of steps
2061 * bDone will be TRUE when the minimization has converged
2062 * bAbort will be TRUE when nsteps steps have been performed or when
2063 * the stepsize becomes smaller than is reasonable for machine precision
2068 while( !bDone && !bAbort ) {
2069 bAbort = (nsteps >= 0) && (count == nsteps);
2071 /* set new coordinates, except for first step */
2073 do_em_step(cr,inputrec,mdatoms,s_min,stepsize,s_min->f,s_try,
2074 constr,top,nrnb,wcycle,count);
2077 evaluate_energy(fplog,bVerbose,cr,
2078 state_global,top_global,s_try,top,
2079 inputrec,nrnb,wcycle,gstat,
2080 vsite,constr,fcd,graph,mdatoms,fr,
2081 mu_tot,enerd,vir,pres,count,count==0);
2084 print_ebin_header(fplog,count,count,s_try->s.lambda[efptFEP]);
2087 s_min->epot = s_try->epot + 1;
2089 /* Print it if necessary */
2092 fprintf(stderr,"Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2093 count,ustep,s_try->epot,s_try->fmax,s_try->a_fmax+1,
2094 (s_try->epot < s_min->epot) ? '\n' : '\r');
2097 if (s_try->epot < s_min->epot) {
2098 /* Store the new (lower) energies */
2099 upd_mdebin(mdebin,FALSE,FALSE,(double)count,
2100 mdatoms->tmass,enerd,&s_try->s,inputrec->fepvals,inputrec->expandedvals,
2101 s_try->s.box, NULL,NULL,vir,pres,NULL,mu_tot,constr);
2102 print_ebin(outf->fp_ene,TRUE,
2103 do_per_step(steps_accepted,inputrec->nstdisreout),
2104 do_per_step(steps_accepted,inputrec->nstorireout),
2105 fplog,count,count,eprNORMAL,TRUE,
2106 mdebin,fcd,&(top_global->groups),&(inputrec->opts));
2111 /* Now if the new energy is smaller than the previous...
2112 * or if this is the first step!
2113 * or if we did random steps!
2116 if ( (count==0) || (s_try->epot < s_min->epot) ) {
2119 /* Test whether the convergence criterion is met... */
2120 bDone = (s_try->fmax < inputrec->em_tol);
2122 /* Copy the arrays for force, positions and energy */
2123 /* The 'Min' array always holds the coords and forces of the minimal
2125 swap_em_state(s_min,s_try);
2129 /* Write to trn, if necessary */
2130 do_x = do_per_step(steps_accepted,inputrec->nstxout);
2131 do_f = do_per_step(steps_accepted,inputrec->nstfout);
2132 write_em_traj(fplog,cr,outf,do_x,do_f,NULL,
2133 top_global,inputrec,count,
2134 s_min,state_global,f_global);
2137 /* If energy is not smaller make the step smaller... */
2140 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count) {
2141 /* Reload the old state */
2142 em_dd_partition_system(fplog,count,cr,top_global,inputrec,
2143 s_min,top,mdatoms,fr,vsite,constr,
2148 /* Determine new step */
2149 stepsize = ustep/s_min->fmax;
2151 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2153 if (count == nsteps || ustep < 1e-12)
2155 if (count == nsteps || ustep < 1e-6)
2160 warn_step(stderr,inputrec->em_tol,count==nsteps,constr!=NULL);
2161 warn_step(fplog ,inputrec->em_tol,count==nsteps,constr!=NULL);
2167 } /* End of the loop */
2169 /* Print some shit... */
2171 fprintf(stderr,"\nwriting lowest energy coordinates.\n");
2172 write_em_traj(fplog,cr,outf,TRUE,inputrec->nstfout,ftp2fn(efSTO,nfile,fnm),
2173 top_global,inputrec,count,
2174 s_min,state_global,f_global);
2176 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2179 print_converged(stderr,SD,inputrec->em_tol,count,bDone,nsteps,
2180 s_min->epot,s_min->fmax,s_min->a_fmax,fnormn);
2181 print_converged(fplog,SD,inputrec->em_tol,count,bDone,nsteps,
2182 s_min->epot,s_min->fmax,s_min->a_fmax,fnormn);
2185 finish_em(fplog,cr,outf,runtime,wcycle);
2187 /* To print the actual number of steps we needed somewhere */
2188 inputrec->nsteps=count;
2190 runtime->nsteps_done = count;
2193 } /* That's all folks */
2196 double do_nm(FILE *fplog,t_commrec *cr,
2197 int nfile,const t_filenm fnm[],
2198 const output_env_t oenv, gmx_bool bVerbose,gmx_bool bCompact,
2200 gmx_vsite_t *vsite,gmx_constr_t constr,
2202 t_inputrec *inputrec,
2203 gmx_mtop_t *top_global,t_fcdata *fcd,
2204 t_state *state_global,
2206 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
2209 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2210 gmx_membed_t membed,
2211 real cpt_period,real max_hours,
2212 const char *deviceOptions,
2213 unsigned long Flags,
2214 gmx_runtime_t *runtime)
2216 const char *NM = "Normal Mode Analysis";
2221 gmx_localtop_t *top;
2222 gmx_enerdata_t *enerd;
2224 gmx_global_stat_t gstat;
2226 real t,t0,lambda,lam0;
2231 gmx_bool bSparse; /* use sparse matrix storage format */
2233 gmx_sparsematrix_t * sparse_matrix = NULL;
2234 real * full_matrix = NULL;
2235 em_state_t * state_work;
2237 /* added with respect to mdrun */
2239 real der_range=10.0*sqrt(GMX_REAL_EPS);
2245 gmx_fatal(FARGS,"Constraints present with Normal Mode Analysis, this combination is not supported");
2248 state_work = init_em_state();
2250 /* Init em and store the local state in state_minimum */
2251 init_em(fplog,NM,cr,inputrec,
2252 state_global,top_global,state_work,&top,
2254 nrnb,mu_tot,fr,&enerd,&graph,mdatoms,&gstat,vsite,constr,
2255 nfile,fnm,&outf,NULL);
2257 natoms = top_global->natoms;
2265 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2266 " which MIGHT not be accurate enough for normal mode analysis.\n"
2267 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2268 " are fairly modest even if you recompile in double precision.\n\n");
2272 /* Check if we can/should use sparse storage format.
2274 * Sparse format is only useful when the Hessian itself is sparse, which it
2275 * will be when we use a cutoff.
2276 * For small systems (n<1000) it is easier to always use full matrix format, though.
2278 if(EEL_FULL(fr->eeltype) || fr->rlist==0.0)
2280 fprintf(stderr,"Non-cutoff electrostatics used, forcing full Hessian format.\n");
2283 else if(top_global->natoms < 1000)
2285 fprintf(stderr,"Small system size (N=%d), using full Hessian format.\n",top_global->natoms);
2290 fprintf(stderr,"Using compressed symmetric sparse Hessian format.\n");
2294 sz = DIM*top_global->natoms;
2296 fprintf(stderr,"Allocating Hessian memory...\n\n");
2300 sparse_matrix=gmx_sparsematrix_init(sz);
2301 sparse_matrix->compressed_symmetric = TRUE;
2305 snew(full_matrix,sz*sz);
2308 /* Initial values */
2309 t0 = inputrec->init_t;
2310 lam0 = inputrec->fepvals->init_lambda;
2318 /* Write start time and temperature */
2319 print_em_start(fplog,cr,runtime,wcycle,NM);
2321 /* fudge nr of steps to nr of atoms */
2322 inputrec->nsteps = natoms*2;
2326 fprintf(stderr,"starting normal mode calculation '%s'\n%d steps.\n\n",
2327 *(top_global->name),(int)inputrec->nsteps);
2330 nnodes = cr->nnodes;
2332 /* Make evaluate_energy do a single node force calculation */
2334 evaluate_energy(fplog,bVerbose,cr,
2335 state_global,top_global,state_work,top,
2336 inputrec,nrnb,wcycle,gstat,
2337 vsite,constr,fcd,graph,mdatoms,fr,
2338 mu_tot,enerd,vir,pres,-1,TRUE);
2339 cr->nnodes = nnodes;
2341 /* if forces are not small, warn user */
2342 get_state_f_norm_max(cr,&(inputrec->opts),mdatoms,state_work);
2346 fprintf(stderr,"Maximum force:%12.5e\n",state_work->fmax);
2347 if (state_work->fmax > 1.0e-3)
2349 fprintf(stderr,"Maximum force probably not small enough to");
2350 fprintf(stderr," ensure that you are in an \nenergy well. ");
2351 fprintf(stderr,"Be aware that negative eigenvalues may occur");
2352 fprintf(stderr," when the\nresulting matrix is diagonalized.\n");
2356 /***********************************************************
2358 * Loop over all pairs in matrix
2360 * do_force called twice. Once with positive and
2361 * once with negative displacement
2363 ************************************************************/
2365 /* Steps are divided one by one over the nodes */
2366 for(atom=cr->nodeid; atom<natoms; atom+=nnodes)
2369 for (d=0; d<DIM; d++)
2371 x_min = state_work->s.x[atom][d];
2373 state_work->s.x[atom][d] = x_min - der_range;
2375 /* Make evaluate_energy do a single node force calculation */
2377 evaluate_energy(fplog,bVerbose,cr,
2378 state_global,top_global,state_work,top,
2379 inputrec,nrnb,wcycle,gstat,
2380 vsite,constr,fcd,graph,mdatoms,fr,
2381 mu_tot,enerd,vir,pres,atom*2,FALSE);
2383 for(i=0; i<natoms; i++)
2385 copy_rvec(state_work->f[i], fneg[i]);
2388 state_work->s.x[atom][d] = x_min + der_range;
2390 evaluate_energy(fplog,bVerbose,cr,
2391 state_global,top_global,state_work,top,
2392 inputrec,nrnb,wcycle,gstat,
2393 vsite,constr,fcd,graph,mdatoms,fr,
2394 mu_tot,enerd,vir,pres,atom*2+1,FALSE);
2395 cr->nnodes = nnodes;
2397 /* x is restored to original */
2398 state_work->s.x[atom][d] = x_min;
2400 for(j=0; j<natoms; j++)
2402 for (k=0; (k<DIM); k++)
2405 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2413 #define mpi_type MPI_DOUBLE
2415 #define mpi_type MPI_FLOAT
2417 MPI_Send(dfdx[0],natoms*DIM,mpi_type,MASTERNODE(cr),cr->nodeid,
2418 cr->mpi_comm_mygroup);
2423 for(node=0; (node<nnodes && atom+node<natoms); node++)
2429 MPI_Recv(dfdx[0],natoms*DIM,mpi_type,node,node,
2430 cr->mpi_comm_mygroup,&stat);
2435 row = (atom + node)*DIM + d;
2437 for(j=0; j<natoms; j++)
2439 for(k=0; k<DIM; k++)
2445 if (col >= row && dfdx[j][k] != 0.0)
2447 gmx_sparsematrix_increment_value(sparse_matrix,
2448 row,col,dfdx[j][k]);
2453 full_matrix[row*sz+col] = dfdx[j][k];
2460 if (bVerbose && fplog)
2465 /* write progress */
2466 if (MASTER(cr) && bVerbose)
2468 fprintf(stderr,"\rFinished step %d out of %d",
2469 min(atom+nnodes,natoms),natoms);
2476 fprintf(stderr,"\n\nWriting Hessian...\n");
2477 gmx_mtxio_write(ftp2fn(efMTX,nfile,fnm),sz,sz,full_matrix,sparse_matrix);
2480 finish_em(fplog,cr,outf,runtime,wcycle);
2482 runtime->nsteps_done = natoms*2;