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-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, 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.
57 #include "gmx_fatal.h"
69 #include "md_support.h"
73 #include "sparsematrix.h"
77 #include "gmx_wallcycle.h"
78 #include "mtop_util.h"
82 #include "gmx_omp_nthreads.h"
83 #include "md_logging.h"
95 static em_state_t *init_em_state()
101 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
102 snew(ems->s.lambda, efptNR);
107 static void print_em_start(FILE *fplog, t_commrec *cr, gmx_runtime_t *runtime,
108 gmx_wallcycle_t wcycle,
113 runtime_start(runtime);
115 sprintf(buf, "Started %s", name);
116 print_date_and_time(fplog, cr->nodeid, buf, NULL);
118 wallcycle_start(wcycle, ewcRUN);
120 static void em_time_end(FILE *fplog, t_commrec *cr, gmx_runtime_t *runtime,
121 gmx_wallcycle_t wcycle)
123 wallcycle_stop(wcycle, ewcRUN);
125 runtime_end(runtime);
128 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
131 fprintf(out, "%s:\n", minimizer);
132 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
133 fprintf(out, " Number of steps = %12d\n", nsteps);
136 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
142 "\nEnergy minimization reached the maximum number "
143 "of steps before the forces reached the requested "
144 "precision Fmax < %g.\n", ftol);
149 "\nEnergy minimization has stopped, but the forces have "
150 "not converged to the requested precision Fmax < %g (which "
151 "may not be possible for your system). It stopped "
152 "because the algorithm tried to make a new step whose size "
153 "was too small, or there was no change in the energy since "
154 "last step. Either way, we regard the minimization as "
155 "converged to within the available machine precision, "
156 "given your starting configuration and EM parameters.\n%s%s",
158 sizeof(real) < sizeof(double) ?
159 "\nDouble precision normally gives you higher accuracy, but "
160 "this is often not needed for preparing to run molecular "
164 "You might need to increase your constraint accuracy, or turn\n"
165 "off constraints altogether (set constraints = none in mdp file)\n" :
168 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
173 static void print_converged(FILE *fp, const char *alg, real ftol,
174 gmx_large_int_t count, gmx_bool bDone, gmx_large_int_t nsteps,
175 real epot, real fmax, int nfmax, real fnorm)
177 char buf[STEPSTRSIZE];
181 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
182 alg, ftol, gmx_step_str(count, buf));
184 else if (count < nsteps)
186 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
187 "but did not reach the requested Fmax < %g.\n",
188 alg, gmx_step_str(count, buf), ftol);
192 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
193 alg, ftol, gmx_step_str(count, buf));
197 fprintf(fp, "Potential Energy = %21.14e\n", epot);
198 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
199 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
201 fprintf(fp, "Potential Energy = %14.7e\n", epot);
202 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
203 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
207 static void get_f_norm_max(t_commrec *cr,
208 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
209 real *fnorm, real *fmax, int *a_fmax)
212 real fmax2, fmax2_0, fam;
213 int la_max, a_max, start, end, i, m, gf;
215 /* This routine finds the largest force and returns it.
216 * On parallel machines the global max is taken.
222 start = mdatoms->start;
223 end = mdatoms->homenr + start;
224 if (mdatoms->cFREEZE)
226 for (i = start; i < end; i++)
228 gf = mdatoms->cFREEZE[i];
230 for (m = 0; m < DIM; m++)
232 if (!opts->nFreeze[gf][m])
247 for (i = start; i < end; i++)
259 if (la_max >= 0 && DOMAINDECOMP(cr))
261 a_max = cr->dd->gatindex[la_max];
269 snew(sum, 2*cr->nnodes+1);
270 sum[2*cr->nodeid] = fmax2;
271 sum[2*cr->nodeid+1] = a_max;
272 sum[2*cr->nnodes] = fnorm2;
273 gmx_sumd(2*cr->nnodes+1, sum, cr);
274 fnorm2 = sum[2*cr->nnodes];
275 /* Determine the global maximum */
276 for (i = 0; i < cr->nnodes; i++)
278 if (sum[2*i] > fmax2)
281 a_max = (int)(sum[2*i+1] + 0.5);
289 *fnorm = sqrt(fnorm2);
301 static void get_state_f_norm_max(t_commrec *cr,
302 t_grpopts *opts, t_mdatoms *mdatoms,
305 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
308 void init_em(FILE *fplog, const char *title,
309 t_commrec *cr, t_inputrec *ir,
310 t_state *state_global, gmx_mtop_t *top_global,
311 em_state_t *ems, gmx_localtop_t **top,
312 rvec **f, rvec **f_global,
313 t_nrnb *nrnb, rvec mu_tot,
314 t_forcerec *fr, gmx_enerdata_t **enerd,
315 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
316 gmx_vsite_t *vsite, gmx_constr_t constr,
317 int nfile, const t_filenm fnm[],
318 gmx_mdoutf_t **outf, t_mdebin **mdebin)
320 int start, homenr, i;
325 fprintf(fplog, "Initiating %s\n", title);
328 state_global->ngtc = 0;
330 /* Initialize lambda variables */
331 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
335 if (DOMAINDECOMP(cr))
337 *top = dd_init_local_top(top_global);
339 dd_init_local_state(cr->dd, state_global, &ems->s);
343 /* Distribute the charge groups over the nodes from the master node */
344 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
345 state_global, top_global, ir,
346 &ems->s, &ems->f, mdatoms, *top,
347 fr, vsite, NULL, constr,
349 dd_store_state(cr->dd, &ems->s);
353 snew(*f_global, top_global->natoms);
363 snew(*f, top_global->natoms);
365 /* Just copy the state */
366 ems->s = *state_global;
367 snew(ems->s.x, ems->s.nalloc);
368 snew(ems->f, ems->s.nalloc);
369 for (i = 0; i < state_global->natoms; i++)
371 copy_rvec(state_global->x[i], ems->s.x[i]);
373 copy_mat(state_global->box, ems->s.box);
375 if (PAR(cr) && ir->eI != eiNM)
377 /* Initialize the particle decomposition and split the topology */
378 *top = split_system(fplog, top_global, ir, cr);
380 pd_cg_range(cr, &fr->cg0, &fr->hcg);
384 *top = gmx_mtop_generate_local_top(top_global, ir);
388 forcerec_set_excl_load(fr, *top, cr);
390 init_bonded_thread_force_reduction(fr, &(*top)->idef);
392 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
394 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
403 pd_at_range(cr, &start, &homenr);
409 homenr = top_global->natoms;
411 atoms2md(top_global, ir, 0, NULL, start, homenr, mdatoms);
412 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
416 set_vsite_top(vsite, *top, mdatoms, cr);
422 if (ir->eConstrAlg == econtSHAKE &&
423 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
425 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
426 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
429 if (!DOMAINDECOMP(cr))
431 set_constraints(constr, *top, ir, mdatoms, cr);
434 if (!ir->bContinuation)
436 /* Constrain the starting coordinates */
438 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
439 ir, NULL, cr, -1, 0, mdatoms,
440 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
441 ems->s.lambda[efptFEP], &dvdl_constr,
442 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
448 *gstat = global_stat_init(ir);
451 *outf = init_mdoutf(nfile, fnm, 0, cr, ir, NULL);
454 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
459 /* Init bin for energy stuff */
460 *mdebin = init_mdebin((*outf)->fp_ene, top_global, ir, NULL);
464 calc_shifts(ems->s.box, fr->shift_vec);
467 static void finish_em(FILE *fplog, t_commrec *cr, gmx_mdoutf_t *outf,
468 gmx_runtime_t *runtime, gmx_wallcycle_t wcycle)
470 if (!(cr->duty & DUTY_PME))
472 /* Tell the PME only node to finish */
473 gmx_pme_send_finish(cr);
478 em_time_end(fplog, cr, runtime, wcycle);
481 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
490 static void copy_em_coords(em_state_t *ems, t_state *state)
494 for (i = 0; (i < state->natoms); i++)
496 copy_rvec(ems->s.x[i], state->x[i]);
500 static void write_em_traj(FILE *fplog, t_commrec *cr,
502 gmx_bool bX, gmx_bool bF, const char *confout,
503 gmx_mtop_t *top_global,
504 t_inputrec *ir, gmx_large_int_t step,
506 t_state *state_global, rvec *f_global)
510 if ((bX || bF || confout != NULL) && !DOMAINDECOMP(cr))
512 copy_em_coords(state, state_global);
519 mdof_flags |= MDOF_X;
523 mdof_flags |= MDOF_F;
525 write_traj(fplog, cr, outf, mdof_flags,
526 top_global, step, (double)step,
527 &state->s, state_global, state->f, f_global, NULL, NULL);
529 if (confout != NULL && MASTER(cr))
531 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
533 /* Make molecules whole only for confout writing */
534 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
538 write_sto_conf_mtop(confout,
539 *top_global->name, top_global,
540 state_global->x, NULL, ir->ePBC, state_global->box);
544 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
546 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
547 gmx_constr_t constr, gmx_localtop_t *top,
548 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
549 gmx_large_int_t count)
561 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
563 gmx_incons("state mismatch in do_em_step");
566 s2->flags = s1->flags;
568 if (s2->nalloc != s1->nalloc)
570 s2->nalloc = s1->nalloc;
571 srenew(s2->x, s1->nalloc);
572 srenew(ems2->f, s1->nalloc);
573 if (s2->flags & (1<<estCGP))
575 srenew(s2->cg_p, s1->nalloc);
579 s2->natoms = s1->natoms;
580 copy_mat(s1->box, s2->box);
581 /* Copy free energy state */
582 for (i = 0; i < efptNR; i++)
584 s2->lambda[i] = s1->lambda[i];
586 copy_mat(s1->box, s2->box);
589 end = md->start + md->homenr;
594 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
599 #pragma omp for schedule(static) nowait
600 for (i = start; i < end; i++)
606 for (m = 0; m < DIM; m++)
608 if (ir->opts.nFreeze[gf][m])
614 x2[i][m] = x1[i][m] + a*f[i][m];
619 if (s2->flags & (1<<estCGP))
621 /* Copy the CG p vector */
624 #pragma omp for schedule(static) nowait
625 for (i = start; i < end; i++)
627 copy_rvec(x1[i], x2[i]);
631 if (DOMAINDECOMP(cr))
633 s2->ddp_count = s1->ddp_count;
634 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
637 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
638 srenew(s2->cg_gl, s2->cg_gl_nalloc);
641 s2->ncg_gl = s1->ncg_gl;
642 #pragma omp for schedule(static) nowait
643 for (i = 0; i < s2->ncg_gl; i++)
645 s2->cg_gl[i] = s1->cg_gl[i];
647 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
653 wallcycle_start(wcycle, ewcCONSTR);
655 constrain(NULL, TRUE, TRUE, constr, &top->idef,
656 ir, NULL, cr, count, 0, md,
657 s1->x, s2->x, NULL, bMolPBC, s2->box,
658 s2->lambda[efptBONDED], &dvdl_constr,
659 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
660 wallcycle_stop(wcycle, ewcCONSTR);
664 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
665 gmx_mtop_t *top_global, t_inputrec *ir,
666 em_state_t *ems, gmx_localtop_t *top,
667 t_mdatoms *mdatoms, t_forcerec *fr,
668 gmx_vsite_t *vsite, gmx_constr_t constr,
669 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
671 /* Repartition the domain decomposition */
672 wallcycle_start(wcycle, ewcDOMDEC);
673 dd_partition_system(fplog, step, cr, FALSE, 1,
674 NULL, top_global, ir,
676 mdatoms, top, fr, vsite, NULL, constr,
677 nrnb, wcycle, FALSE);
678 dd_store_state(cr->dd, &ems->s);
679 wallcycle_stop(wcycle, ewcDOMDEC);
682 static void evaluate_energy(FILE *fplog, gmx_bool bVerbose, t_commrec *cr,
683 t_state *state_global, gmx_mtop_t *top_global,
684 em_state_t *ems, gmx_localtop_t *top,
685 t_inputrec *inputrec,
686 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
687 gmx_global_stat_t gstat,
688 gmx_vsite_t *vsite, gmx_constr_t constr,
690 t_graph *graph, t_mdatoms *mdatoms,
691 t_forcerec *fr, rvec mu_tot,
692 gmx_enerdata_t *enerd, tensor vir, tensor pres,
693 gmx_large_int_t count, gmx_bool bFirst)
698 tensor force_vir, shake_vir, ekin;
699 real dvdl_constr, prescorr, enercorr, dvdlcorr;
702 /* Set the time to the initial time, the time does not change during EM */
703 t = inputrec->init_t;
706 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
708 /* This the first state or an old state used before the last ns */
714 if (inputrec->nstlist > 0)
718 else if (inputrec->nstlist == -1)
720 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
723 gmx_sumi(1, &nabnsb, cr);
731 construct_vsites(fplog, vsite, ems->s.x, nrnb, 1, NULL,
732 top->idef.iparams, top->idef.il,
733 fr->ePBC, fr->bMolPBC, graph, cr, ems->s.box);
736 if (DOMAINDECOMP(cr))
740 /* Repartition the domain decomposition */
741 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
742 ems, top, mdatoms, fr, vsite, constr,
747 /* Calc force & energy on new trial position */
748 /* do_force always puts the charge groups in the box and shifts again
749 * We do not unshift, so molecules are always whole in congrad.c
751 do_force(fplog, cr, inputrec,
752 count, nrnb, wcycle, top, top_global, &top_global->groups,
753 ems->s.box, ems->s.x, &ems->s.hist,
754 ems->f, force_vir, mdatoms, enerd, fcd,
755 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
756 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
757 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
758 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
760 /* Clear the unused shake virial and pressure */
761 clear_mat(shake_vir);
764 /* Communicate stuff when parallel */
765 if (PAR(cr) && inputrec->eI != eiNM)
767 wallcycle_start(wcycle, ewcMoveE);
769 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
770 inputrec, NULL, NULL, NULL, 1, &terminate,
771 top_global, &ems->s, FALSE,
777 wallcycle_stop(wcycle, ewcMoveE);
780 /* Calculate long range corrections to pressure and energy */
781 calc_dispcorr(fplog, inputrec, fr, count, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
782 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
783 enerd->term[F_DISPCORR] = enercorr;
784 enerd->term[F_EPOT] += enercorr;
785 enerd->term[F_PRES] += prescorr;
786 enerd->term[F_DVDL] += dvdlcorr;
788 ems->epot = enerd->term[F_EPOT];
792 /* Project out the constraint components of the force */
793 wallcycle_start(wcycle, ewcCONSTR);
795 constrain(NULL, FALSE, FALSE, constr, &top->idef,
796 inputrec, NULL, cr, count, 0, mdatoms,
797 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
798 ems->s.lambda[efptBONDED], &dvdl_constr,
799 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
800 if (fr->bSepDVDL && fplog)
802 fprintf(fplog, sepdvdlformat, "Constraints", t, dvdl_constr);
804 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
805 m_add(force_vir, shake_vir, vir);
806 wallcycle_stop(wcycle, ewcCONSTR);
810 copy_mat(force_vir, vir);
814 enerd->term[F_PRES] =
815 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
817 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
819 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
821 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
825 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
827 em_state_t *s_min, em_state_t *s_b)
831 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
833 unsigned char *grpnrFREEZE;
837 fprintf(debug, "Doing reorder_partsum\n");
843 cgs_gl = dd_charge_groups_global(cr->dd);
844 index = cgs_gl->index;
846 /* Collect fm in a global vector fmg.
847 * This conflicts with the spirit of domain decomposition,
848 * but to fully optimize this a much more complicated algorithm is required.
850 snew(fmg, mtop->natoms);
852 ncg = s_min->s.ncg_gl;
853 cg_gl = s_min->s.cg_gl;
855 for (c = 0; c < ncg; c++)
860 for (a = a0; a < a1; a++)
862 copy_rvec(fm[i], fmg[a]);
866 gmx_sum(mtop->natoms*3, fmg[0], cr);
868 /* Now we will determine the part of the sum for the cgs in state s_b */
870 cg_gl = s_b->s.cg_gl;
874 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
875 for (c = 0; c < ncg; c++)
880 for (a = a0; a < a1; a++)
882 if (mdatoms->cFREEZE && grpnrFREEZE)
886 for (m = 0; m < DIM; m++)
888 if (!opts->nFreeze[gf][m])
890 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
902 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
904 em_state_t *s_min, em_state_t *s_b)
910 /* This is just the classical Polak-Ribiere calculation of beta;
911 * it looks a bit complicated since we take freeze groups into account,
912 * and might have to sum it in parallel runs.
915 if (!DOMAINDECOMP(cr) ||
916 (s_min->s.ddp_count == cr->dd->ddp_count &&
917 s_b->s.ddp_count == cr->dd->ddp_count))
923 /* This part of code can be incorrect with DD,
924 * since the atom ordering in s_b and s_min might differ.
926 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
928 if (mdatoms->cFREEZE)
930 gf = mdatoms->cFREEZE[i];
932 for (m = 0; m < DIM; m++)
934 if (!opts->nFreeze[gf][m])
936 sum += (fb[i][m] - fm[i][m])*fb[i][m];
943 /* We need to reorder cgs while summing */
944 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
948 gmx_sumd(1, &sum, cr);
951 return sum/sqr(s_min->fnorm);
954 double do_cg(FILE *fplog, t_commrec *cr,
955 int nfile, const t_filenm fnm[],
956 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
958 gmx_vsite_t *vsite, gmx_constr_t constr,
960 t_inputrec *inputrec,
961 gmx_mtop_t *top_global, t_fcdata *fcd,
962 t_state *state_global,
964 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
967 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
969 real cpt_period, real max_hours,
970 const char *deviceOptions,
972 gmx_runtime_t *runtime)
974 const char *CG = "Polak-Ribiere Conjugate Gradients";
976 em_state_t *s_min, *s_a, *s_b, *s_c;
978 gmx_enerdata_t *enerd;
980 gmx_global_stat_t gstat;
982 rvec *f_global, *p, *sf, *sfm;
983 double gpa, gpb, gpc, tmp, sum[2], minstep;
986 real a, b, c, beta = 0.0;
990 gmx_bool converged, foundlower;
992 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
994 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
996 int i, m, gf, step, nminstep;
1001 s_min = init_em_state();
1002 s_a = init_em_state();
1003 s_b = init_em_state();
1004 s_c = init_em_state();
1006 /* Init em and store the local state in s_min */
1007 init_em(fplog, CG, cr, inputrec,
1008 state_global, top_global, s_min, &top, &f, &f_global,
1009 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1010 nfile, fnm, &outf, &mdebin);
1012 /* Print to log file */
1013 print_em_start(fplog, cr, runtime, wcycle, CG);
1015 /* Max number of steps */
1016 number_steps = inputrec->nsteps;
1020 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1024 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1027 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1028 /* do_force always puts the charge groups in the box and shifts again
1029 * We do not unshift, so molecules are always whole in congrad.c
1031 evaluate_energy(fplog, bVerbose, cr,
1032 state_global, top_global, s_min, top,
1033 inputrec, nrnb, wcycle, gstat,
1034 vsite, constr, fcd, graph, mdatoms, fr,
1035 mu_tot, enerd, vir, pres, -1, TRUE);
1040 /* Copy stuff to the energy bin for easy printing etc. */
1041 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1042 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1043 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1045 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1046 print_ebin(outf->fp_ene, TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1047 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1051 /* Estimate/guess the initial stepsize */
1052 stepsize = inputrec->em_stepsize/s_min->fnorm;
1056 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1057 s_min->fmax, s_min->a_fmax+1);
1058 fprintf(stderr, " F-Norm = %12.5e\n",
1059 s_min->fnorm/sqrt(state_global->natoms));
1060 fprintf(stderr, "\n");
1061 /* and copy to the log file too... */
1062 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1063 s_min->fmax, s_min->a_fmax+1);
1064 fprintf(fplog, " F-Norm = %12.5e\n",
1065 s_min->fnorm/sqrt(state_global->natoms));
1066 fprintf(fplog, "\n");
1068 /* Start the loop over CG steps.
1069 * Each successful step is counted, and we continue until
1070 * we either converge or reach the max number of steps.
1073 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1076 /* start taking steps in a new direction
1077 * First time we enter the routine, beta=0, and the direction is
1078 * simply the negative gradient.
1081 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1086 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1088 if (mdatoms->cFREEZE)
1090 gf = mdatoms->cFREEZE[i];
1092 for (m = 0; m < DIM; m++)
1094 if (!inputrec->opts.nFreeze[gf][m])
1096 p[i][m] = sf[i][m] + beta*p[i][m];
1097 gpa -= p[i][m]*sf[i][m];
1098 /* f is negative gradient, thus the sign */
1107 /* Sum the gradient along the line across CPUs */
1110 gmx_sumd(1, &gpa, cr);
1113 /* Calculate the norm of the search vector */
1114 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1116 /* Just in case stepsize reaches zero due to numerical precision... */
1119 stepsize = inputrec->em_stepsize/pnorm;
1123 * Double check the value of the derivative in the search direction.
1124 * If it is positive it must be due to the old information in the
1125 * CG formula, so just remove that and start over with beta=0.
1126 * This corresponds to a steepest descent step.
1131 step--; /* Don't count this step since we are restarting */
1132 continue; /* Go back to the beginning of the big for-loop */
1135 /* Calculate minimum allowed stepsize, before the average (norm)
1136 * relative change in coordinate is smaller than precision
1139 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1141 for (m = 0; m < DIM; m++)
1143 tmp = fabs(s_min->s.x[i][m]);
1152 /* Add up from all CPUs */
1155 gmx_sumd(1, &minstep, cr);
1158 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1160 if (stepsize < minstep)
1166 /* Write coordinates if necessary */
1167 do_x = do_per_step(step, inputrec->nstxout);
1168 do_f = do_per_step(step, inputrec->nstfout);
1170 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1171 top_global, inputrec, step,
1172 s_min, state_global, f_global);
1174 /* Take a step downhill.
1175 * In theory, we should minimize the function along this direction.
1176 * That is quite possible, but it turns out to take 5-10 function evaluations
1177 * for each line. However, we dont really need to find the exact minimum -
1178 * it is much better to start a new CG step in a modified direction as soon
1179 * as we are close to it. This will save a lot of energy evaluations.
1181 * In practice, we just try to take a single step.
1182 * If it worked (i.e. lowered the energy), we increase the stepsize but
1183 * the continue straight to the next CG step without trying to find any minimum.
1184 * If it didn't work (higher energy), there must be a minimum somewhere between
1185 * the old position and the new one.
1187 * Due to the finite numerical accuracy, it turns out that it is a good idea
1188 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1189 * This leads to lower final energies in the tests I've done. / Erik
1191 s_a->epot = s_min->epot;
1193 c = a + stepsize; /* reference position along line is zero */
1195 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1197 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1198 s_min, top, mdatoms, fr, vsite, constr,
1202 /* Take a trial step (new coords in s_c) */
1203 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1204 constr, top, nrnb, wcycle, -1);
1207 /* Calculate energy for the trial step */
1208 evaluate_energy(fplog, bVerbose, cr,
1209 state_global, top_global, s_c, top,
1210 inputrec, nrnb, wcycle, gstat,
1211 vsite, constr, fcd, graph, mdatoms, fr,
1212 mu_tot, enerd, vir, pres, -1, FALSE);
1214 /* Calc derivative along line */
1218 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1220 for (m = 0; m < DIM; m++)
1222 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1225 /* Sum the gradient along the line across CPUs */
1228 gmx_sumd(1, &gpc, cr);
1231 /* This is the max amount of increase in energy we tolerate */
1232 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1234 /* Accept the step if the energy is lower, or if it is not significantly higher
1235 * and the line derivative is still negative.
1237 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1240 /* Great, we found a better energy. Increase step for next iteration
1241 * if we are still going down, decrease it otherwise
1245 stepsize *= 1.618034; /* The golden section */
1249 stepsize *= 0.618034; /* 1/golden section */
1254 /* New energy is the same or higher. We will have to do some work
1255 * to find a smaller value in the interval. Take smaller step next time!
1258 stepsize *= 0.618034;
1264 /* OK, if we didn't find a lower value we will have to locate one now - there must
1265 * be one in the interval [a=0,c].
1266 * The same thing is valid here, though: Don't spend dozens of iterations to find
1267 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1268 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1270 * I also have a safeguard for potentially really patological functions so we never
1271 * take more than 20 steps before we give up ...
1273 * If we already found a lower value we just skip this step and continue to the update.
1281 /* Select a new trial point.
1282 * If the derivatives at points a & c have different sign we interpolate to zero,
1283 * otherwise just do a bisection.
1285 if (gpa < 0 && gpc > 0)
1287 b = a + gpa*(a-c)/(gpc-gpa);
1294 /* safeguard if interpolation close to machine accuracy causes errors:
1295 * never go outside the interval
1297 if (b <= a || b >= c)
1302 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1304 /* Reload the old state */
1305 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1306 s_min, top, mdatoms, fr, vsite, constr,
1310 /* Take a trial step to this new point - new coords in s_b */
1311 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1312 constr, top, nrnb, wcycle, -1);
1315 /* Calculate energy for the trial step */
1316 evaluate_energy(fplog, bVerbose, cr,
1317 state_global, top_global, s_b, top,
1318 inputrec, nrnb, wcycle, gstat,
1319 vsite, constr, fcd, graph, mdatoms, fr,
1320 mu_tot, enerd, vir, pres, -1, FALSE);
1322 /* p does not change within a step, but since the domain decomposition
1323 * might change, we have to use cg_p of s_b here.
1328 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1330 for (m = 0; m < DIM; m++)
1332 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1335 /* Sum the gradient along the line across CPUs */
1338 gmx_sumd(1, &gpb, cr);
1343 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1344 s_a->epot, s_b->epot, s_c->epot, gpb);
1347 epot_repl = s_b->epot;
1349 /* Keep one of the intervals based on the value of the derivative at the new point */
1352 /* Replace c endpoint with b */
1353 swap_em_state(s_b, s_c);
1359 /* Replace a endpoint with b */
1360 swap_em_state(s_b, s_a);
1366 * Stop search as soon as we find a value smaller than the endpoints.
1367 * Never run more than 20 steps, no matter what.
1371 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1374 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1377 /* OK. We couldn't find a significantly lower energy.
1378 * If beta==0 this was steepest descent, and then we give up.
1379 * If not, set beta=0 and restart with steepest descent before quitting.
1389 /* Reset memory before giving up */
1395 /* Select min energy state of A & C, put the best in B.
1397 if (s_c->epot < s_a->epot)
1401 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1402 s_c->epot, s_a->epot);
1404 swap_em_state(s_b, s_c);
1412 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1413 s_a->epot, s_c->epot);
1415 swap_em_state(s_b, s_a);
1425 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1428 swap_em_state(s_b, s_c);
1433 /* new search direction */
1434 /* beta = 0 means forget all memory and restart with steepest descents. */
1435 if (nstcg && ((step % nstcg) == 0))
1441 /* s_min->fnorm cannot be zero, because then we would have converged
1445 /* Polak-Ribiere update.
1446 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1448 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1450 /* Limit beta to prevent oscillations */
1451 if (fabs(beta) > 5.0)
1457 /* update positions */
1458 swap_em_state(s_min, s_b);
1461 /* Print it if necessary */
1466 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1467 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1468 s_min->fmax, s_min->a_fmax+1);
1470 /* Store the new (lower) energies */
1471 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1472 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1473 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1475 do_log = do_per_step(step, inputrec->nstlog);
1476 do_ene = do_per_step(step, inputrec->nstenergy);
1479 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1481 print_ebin(outf->fp_ene, do_ene, FALSE, FALSE,
1482 do_log ? fplog : NULL, step, step, eprNORMAL,
1483 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1486 /* Stop when the maximum force lies below tolerance.
1487 * If we have reached machine precision, converged is already set to true.
1489 converged = converged || (s_min->fmax < inputrec->em_tol);
1491 } /* End of the loop */
1495 step--; /* we never took that last step in this case */
1498 if (s_min->fmax > inputrec->em_tol)
1502 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1503 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1510 /* If we printed energy and/or logfile last step (which was the last step)
1511 * we don't have to do it again, but otherwise print the final values.
1515 /* Write final value to log since we didn't do anything the last step */
1516 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1518 if (!do_ene || !do_log)
1520 /* Write final energy file entries */
1521 print_ebin(outf->fp_ene, !do_ene, FALSE, FALSE,
1522 !do_log ? fplog : NULL, step, step, eprNORMAL,
1523 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1527 /* Print some stuff... */
1530 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1534 * For accurate normal mode calculation it is imperative that we
1535 * store the last conformation into the full precision binary trajectory.
1537 * However, we should only do it if we did NOT already write this step
1538 * above (which we did if do_x or do_f was true).
1540 do_x = !do_per_step(step, inputrec->nstxout);
1541 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1543 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1544 top_global, inputrec, step,
1545 s_min, state_global, f_global);
1547 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1551 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1552 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1553 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1554 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1556 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1559 finish_em(fplog, cr, outf, runtime, wcycle);
1561 /* To print the actual number of steps we needed somewhere */
1562 runtime->nsteps_done = step;
1565 } /* That's all folks */
1568 double do_lbfgs(FILE *fplog, t_commrec *cr,
1569 int nfile, const t_filenm fnm[],
1570 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
1572 gmx_vsite_t *vsite, gmx_constr_t constr,
1574 t_inputrec *inputrec,
1575 gmx_mtop_t *top_global, t_fcdata *fcd,
1578 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1581 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1582 gmx_membed_t membed,
1583 real cpt_period, real max_hours,
1584 const char *deviceOptions,
1585 unsigned long Flags,
1586 gmx_runtime_t *runtime)
1588 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1590 gmx_localtop_t *top;
1591 gmx_enerdata_t *enerd;
1593 gmx_global_stat_t gstat;
1596 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1597 double stepsize, gpa, gpb, gpc, tmp, minstep;
1598 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1599 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1600 real a, b, c, maxdelta, delta;
1601 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1602 real dgdx, dgdg, sq, yr, beta;
1604 gmx_bool converged, first;
1607 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1609 int start, end, number_steps;
1611 int i, k, m, n, nfmax, gf, step;
1618 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1623 gmx_fatal(FARGS, "The combination of constraints and L-BFGS minimization is not implemented. Either do not use constraints, or use another minimizer (e.g. steepest descent).");
1626 n = 3*state->natoms;
1627 nmaxcorr = inputrec->nbfgscorr;
1629 /* Allocate memory */
1630 /* Use pointers to real so we dont have to loop over both atoms and
1631 * dimensions all the time...
1632 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1633 * that point to the same memory.
1646 snew(rho, nmaxcorr);
1647 snew(alpha, nmaxcorr);
1650 for (i = 0; i < nmaxcorr; i++)
1656 for (i = 0; i < nmaxcorr; i++)
1665 init_em(fplog, LBFGS, cr, inputrec,
1666 state, top_global, &ems, &top, &f, &f_global,
1667 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1668 nfile, fnm, &outf, &mdebin);
1669 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1670 * so we free some memory again.
1675 xx = (real *)state->x;
1678 start = mdatoms->start;
1679 end = mdatoms->homenr + start;
1681 /* Print to log file */
1682 print_em_start(fplog, cr, runtime, wcycle, LBFGS);
1684 do_log = do_ene = do_x = do_f = TRUE;
1686 /* Max number of steps */
1687 number_steps = inputrec->nsteps;
1689 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1691 for (i = start; i < end; i++)
1693 if (mdatoms->cFREEZE)
1695 gf = mdatoms->cFREEZE[i];
1697 for (m = 0; m < DIM; m++)
1699 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1704 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1708 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1713 construct_vsites(fplog, vsite, state->x, nrnb, 1, NULL,
1714 top->idef.iparams, top->idef.il,
1715 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1718 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1719 /* do_force always puts the charge groups in the box and shifts again
1720 * We do not unshift, so molecules are always whole
1725 evaluate_energy(fplog, bVerbose, cr,
1726 state, top_global, &ems, top,
1727 inputrec, nrnb, wcycle, gstat,
1728 vsite, constr, fcd, graph, mdatoms, fr,
1729 mu_tot, enerd, vir, pres, -1, TRUE);
1734 /* Copy stuff to the energy bin for easy printing etc. */
1735 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1736 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1737 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1739 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1740 print_ebin(outf->fp_ene, TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1741 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1745 /* This is the starting energy */
1746 Epot = enerd->term[F_EPOT];
1752 /* Set the initial step.
1753 * since it will be multiplied by the non-normalized search direction
1754 * vector (force vector the first time), we scale it by the
1755 * norm of the force.
1760 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1761 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1762 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1763 fprintf(stderr, "\n");
1764 /* and copy to the log file too... */
1765 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1766 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1767 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1768 fprintf(fplog, "\n");
1772 for (i = 0; i < n; i++)
1776 dx[point][i] = ff[i]; /* Initial search direction */
1784 stepsize = 1.0/fnorm;
1787 /* Start the loop over BFGS steps.
1788 * Each successful step is counted, and we continue until
1789 * we either converge or reach the max number of steps.
1794 /* Set the gradient from the force */
1796 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1799 /* Write coordinates if necessary */
1800 do_x = do_per_step(step, inputrec->nstxout);
1801 do_f = do_per_step(step, inputrec->nstfout);
1806 mdof_flags |= MDOF_X;
1811 mdof_flags |= MDOF_F;
1814 write_traj(fplog, cr, outf, mdof_flags,
1815 top_global, step, (real)step, state, state, f, f, NULL, NULL);
1817 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1819 /* pointer to current direction - point=0 first time here */
1822 /* calculate line gradient */
1823 for (gpa = 0, i = 0; i < n; i++)
1828 /* Calculate minimum allowed stepsize, before the average (norm)
1829 * relative change in coordinate is smaller than precision
1831 for (minstep = 0, i = 0; i < n; i++)
1841 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1843 if (stepsize < minstep)
1849 /* Store old forces and coordinates */
1850 for (i = 0; i < n; i++)
1859 for (i = 0; i < n; i++)
1864 /* Take a step downhill.
1865 * In theory, we should minimize the function along this direction.
1866 * That is quite possible, but it turns out to take 5-10 function evaluations
1867 * for each line. However, we dont really need to find the exact minimum -
1868 * it is much better to start a new BFGS step in a modified direction as soon
1869 * as we are close to it. This will save a lot of energy evaluations.
1871 * In practice, we just try to take a single step.
1872 * If it worked (i.e. lowered the energy), we increase the stepsize but
1873 * the continue straight to the next BFGS step without trying to find any minimum.
1874 * If it didn't work (higher energy), there must be a minimum somewhere between
1875 * the old position and the new one.
1877 * Due to the finite numerical accuracy, it turns out that it is a good idea
1878 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1879 * This leads to lower final energies in the tests I've done. / Erik
1884 c = a + stepsize; /* reference position along line is zero */
1886 /* Check stepsize first. We do not allow displacements
1887 * larger than emstep.
1893 for (i = 0; i < n; i++)
1896 if (delta > maxdelta)
1901 if (maxdelta > inputrec->em_stepsize)
1906 while (maxdelta > inputrec->em_stepsize);
1908 /* Take a trial step */
1909 for (i = 0; i < n; i++)
1911 xc[i] = lastx[i] + c*s[i];
1915 /* Calculate energy for the trial step */
1916 ems.s.x = (rvec *)xc;
1918 evaluate_energy(fplog, bVerbose, cr,
1919 state, top_global, &ems, top,
1920 inputrec, nrnb, wcycle, gstat,
1921 vsite, constr, fcd, graph, mdatoms, fr,
1922 mu_tot, enerd, vir, pres, step, FALSE);
1925 /* Calc derivative along line */
1926 for (gpc = 0, i = 0; i < n; i++)
1928 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1930 /* Sum the gradient along the line across CPUs */
1933 gmx_sumd(1, &gpc, cr);
1936 /* This is the max amount of increase in energy we tolerate */
1937 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1939 /* Accept the step if the energy is lower, or if it is not significantly higher
1940 * and the line derivative is still negative.
1942 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1945 /* Great, we found a better energy. Increase step for next iteration
1946 * if we are still going down, decrease it otherwise
1950 stepsize *= 1.618034; /* The golden section */
1954 stepsize *= 0.618034; /* 1/golden section */
1959 /* New energy is the same or higher. We will have to do some work
1960 * to find a smaller value in the interval. Take smaller step next time!
1963 stepsize *= 0.618034;
1966 /* OK, if we didn't find a lower value we will have to locate one now - there must
1967 * be one in the interval [a=0,c].
1968 * The same thing is valid here, though: Don't spend dozens of iterations to find
1969 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1970 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1972 * I also have a safeguard for potentially really patological functions so we never
1973 * take more than 20 steps before we give up ...
1975 * If we already found a lower value we just skip this step and continue to the update.
1984 /* Select a new trial point.
1985 * If the derivatives at points a & c have different sign we interpolate to zero,
1986 * otherwise just do a bisection.
1989 if (gpa < 0 && gpc > 0)
1991 b = a + gpa*(a-c)/(gpc-gpa);
1998 /* safeguard if interpolation close to machine accuracy causes errors:
1999 * never go outside the interval
2001 if (b <= a || b >= c)
2006 /* Take a trial step */
2007 for (i = 0; i < n; i++)
2009 xb[i] = lastx[i] + b*s[i];
2013 /* Calculate energy for the trial step */
2014 ems.s.x = (rvec *)xb;
2016 evaluate_energy(fplog, bVerbose, cr,
2017 state, top_global, &ems, top,
2018 inputrec, nrnb, wcycle, gstat,
2019 vsite, constr, fcd, graph, mdatoms, fr,
2020 mu_tot, enerd, vir, pres, step, FALSE);
2025 for (gpb = 0, i = 0; i < n; i++)
2027 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2030 /* Sum the gradient along the line across CPUs */
2033 gmx_sumd(1, &gpb, cr);
2036 /* Keep one of the intervals based on the value of the derivative at the new point */
2039 /* Replace c endpoint with b */
2043 /* swap coord pointers b/c */
2053 /* Replace a endpoint with b */
2057 /* swap coord pointers a/b */
2067 * Stop search as soon as we find a value smaller than the endpoints,
2068 * or if the tolerance is below machine precision.
2069 * Never run more than 20 steps, no matter what.
2073 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2075 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2077 /* OK. We couldn't find a significantly lower energy.
2078 * If ncorr==0 this was steepest descent, and then we give up.
2079 * If not, reset memory to restart as steepest descent before quitting.
2091 /* Search in gradient direction */
2092 for (i = 0; i < n; i++)
2094 dx[point][i] = ff[i];
2096 /* Reset stepsize */
2097 stepsize = 1.0/fnorm;
2102 /* Select min energy state of A & C, put the best in xx/ff/Epot
2108 for (i = 0; i < n; i++)
2119 for (i = 0; i < n; i++)
2133 for (i = 0; i < n; i++)
2141 /* Update the memory information, and calculate a new
2142 * approximation of the inverse hessian
2145 /* Have new data in Epot, xx, ff */
2146 if (ncorr < nmaxcorr)
2151 for (i = 0; i < n; i++)
2153 dg[point][i] = lastf[i]-ff[i];
2154 dx[point][i] *= stepsize;
2159 for (i = 0; i < n; i++)
2161 dgdg += dg[point][i]*dg[point][i];
2162 dgdx += dg[point][i]*dx[point][i];
2167 rho[point] = 1.0/dgdx;
2170 if (point >= nmaxcorr)
2176 for (i = 0; i < n; i++)
2183 /* Recursive update. First go back over the memory points */
2184 for (k = 0; k < ncorr; k++)
2193 for (i = 0; i < n; i++)
2195 sq += dx[cp][i]*p[i];
2198 alpha[cp] = rho[cp]*sq;
2200 for (i = 0; i < n; i++)
2202 p[i] -= alpha[cp]*dg[cp][i];
2206 for (i = 0; i < n; i++)
2211 /* And then go forward again */
2212 for (k = 0; k < ncorr; k++)
2215 for (i = 0; i < n; i++)
2217 yr += p[i]*dg[cp][i];
2221 beta = alpha[cp]-beta;
2223 for (i = 0; i < n; i++)
2225 p[i] += beta*dx[cp][i];
2235 for (i = 0; i < n; i++)
2239 dx[point][i] = p[i];
2249 /* Test whether the convergence criterion is met */
2250 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2252 /* Print it if necessary */
2257 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2258 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2260 /* Store the new (lower) energies */
2261 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2262 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2263 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2264 do_log = do_per_step(step, inputrec->nstlog);
2265 do_ene = do_per_step(step, inputrec->nstenergy);
2268 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2270 print_ebin(outf->fp_ene, do_ene, FALSE, FALSE,
2271 do_log ? fplog : NULL, step, step, eprNORMAL,
2272 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2275 /* Stop when the maximum force lies below tolerance.
2276 * If we have reached machine precision, converged is already set to true.
2279 converged = converged || (fmax < inputrec->em_tol);
2281 } /* End of the loop */
2285 step--; /* we never took that last step in this case */
2288 if (fmax > inputrec->em_tol)
2292 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2293 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2298 /* If we printed energy and/or logfile last step (which was the last step)
2299 * we don't have to do it again, but otherwise print the final values.
2301 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2303 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2305 if (!do_ene || !do_log) /* Write final energy file entries */
2307 print_ebin(outf->fp_ene, !do_ene, FALSE, FALSE,
2308 !do_log ? fplog : NULL, step, step, eprNORMAL,
2309 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2312 /* Print some stuff... */
2315 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2319 * For accurate normal mode calculation it is imperative that we
2320 * store the last conformation into the full precision binary trajectory.
2322 * However, we should only do it if we did NOT already write this step
2323 * above (which we did if do_x or do_f was true).
2325 do_x = !do_per_step(step, inputrec->nstxout);
2326 do_f = !do_per_step(step, inputrec->nstfout);
2327 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2328 top_global, inputrec, step,
2333 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2334 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2335 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2336 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2338 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2341 finish_em(fplog, cr, outf, runtime, wcycle);
2343 /* To print the actual number of steps we needed somewhere */
2344 runtime->nsteps_done = step;
2347 } /* That's all folks */
2350 double do_steep(FILE *fplog, t_commrec *cr,
2351 int nfile, const t_filenm fnm[],
2352 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
2354 gmx_vsite_t *vsite, gmx_constr_t constr,
2356 t_inputrec *inputrec,
2357 gmx_mtop_t *top_global, t_fcdata *fcd,
2358 t_state *state_global,
2360 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2363 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2364 gmx_membed_t membed,
2365 real cpt_period, real max_hours,
2366 const char *deviceOptions,
2367 unsigned long Flags,
2368 gmx_runtime_t *runtime)
2370 const char *SD = "Steepest Descents";
2371 em_state_t *s_min, *s_try;
2373 gmx_localtop_t *top;
2374 gmx_enerdata_t *enerd;
2376 gmx_global_stat_t gstat;
2378 real stepsize, constepsize;
2382 gmx_bool bDone, bAbort, do_x, do_f;
2387 int steps_accepted = 0;
2391 s_min = init_em_state();
2392 s_try = init_em_state();
2394 /* Init em and store the local state in s_try */
2395 init_em(fplog, SD, cr, inputrec,
2396 state_global, top_global, s_try, &top, &f, &f_global,
2397 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2398 nfile, fnm, &outf, &mdebin);
2400 /* Print to log file */
2401 print_em_start(fplog, cr, runtime, wcycle, SD);
2403 /* Set variables for stepsize (in nm). This is the largest
2404 * step that we are going to make in any direction.
2406 ustep = inputrec->em_stepsize;
2409 /* Max number of steps */
2410 nsteps = inputrec->nsteps;
2414 /* Print to the screen */
2415 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2419 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2422 /**** HERE STARTS THE LOOP ****
2423 * count is the counter for the number of steps
2424 * bDone will be TRUE when the minimization has converged
2425 * bAbort will be TRUE when nsteps steps have been performed or when
2426 * the stepsize becomes smaller than is reasonable for machine precision
2431 while (!bDone && !bAbort)
2433 bAbort = (nsteps >= 0) && (count == nsteps);
2435 /* set new coordinates, except for first step */
2438 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2439 s_min, stepsize, s_min->f, s_try,
2440 constr, top, nrnb, wcycle, count);
2443 evaluate_energy(fplog, bVerbose, cr,
2444 state_global, top_global, s_try, top,
2445 inputrec, nrnb, wcycle, gstat,
2446 vsite, constr, fcd, graph, mdatoms, fr,
2447 mu_tot, enerd, vir, pres, count, count == 0);
2451 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2456 s_min->epot = s_try->epot + 1;
2459 /* Print it if necessary */
2464 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2465 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2466 (s_try->epot < s_min->epot) ? '\n' : '\r');
2469 if (s_try->epot < s_min->epot)
2471 /* Store the new (lower) energies */
2472 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2473 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2474 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2475 print_ebin(outf->fp_ene, TRUE,
2476 do_per_step(steps_accepted, inputrec->nstdisreout),
2477 do_per_step(steps_accepted, inputrec->nstorireout),
2478 fplog, count, count, eprNORMAL, TRUE,
2479 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2484 /* Now if the new energy is smaller than the previous...
2485 * or if this is the first step!
2486 * or if we did random steps!
2489 if ( (count == 0) || (s_try->epot < s_min->epot) )
2493 /* Test whether the convergence criterion is met... */
2494 bDone = (s_try->fmax < inputrec->em_tol);
2496 /* Copy the arrays for force, positions and energy */
2497 /* The 'Min' array always holds the coords and forces of the minimal
2499 swap_em_state(s_min, s_try);
2505 /* Write to trn, if necessary */
2506 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2507 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2508 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2509 top_global, inputrec, count,
2510 s_min, state_global, f_global);
2514 /* If energy is not smaller make the step smaller... */
2517 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2519 /* Reload the old state */
2520 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2521 s_min, top, mdatoms, fr, vsite, constr,
2526 /* Determine new step */
2527 stepsize = ustep/s_min->fmax;
2529 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2531 if (count == nsteps || ustep < 1e-12)
2533 if (count == nsteps || ustep < 1e-6)
2538 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2539 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2545 } /* End of the loop */
2547 /* Print some shit... */
2550 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2552 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2553 top_global, inputrec, count,
2554 s_min, state_global, f_global);
2556 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2560 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2561 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2562 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2563 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2566 finish_em(fplog, cr, outf, runtime, wcycle);
2568 /* To print the actual number of steps we needed somewhere */
2569 inputrec->nsteps = count;
2571 runtime->nsteps_done = count;
2574 } /* That's all folks */
2577 double do_nm(FILE *fplog, t_commrec *cr,
2578 int nfile, const t_filenm fnm[],
2579 const output_env_t oenv, gmx_bool bVerbose, gmx_bool bCompact,
2581 gmx_vsite_t *vsite, gmx_constr_t constr,
2583 t_inputrec *inputrec,
2584 gmx_mtop_t *top_global, t_fcdata *fcd,
2585 t_state *state_global,
2587 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2590 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2591 gmx_membed_t membed,
2592 real cpt_period, real max_hours,
2593 const char *deviceOptions,
2594 unsigned long Flags,
2595 gmx_runtime_t *runtime)
2597 const char *NM = "Normal Mode Analysis";
2599 int natoms, atom, d;
2602 gmx_localtop_t *top;
2603 gmx_enerdata_t *enerd;
2605 gmx_global_stat_t gstat;
2607 real t, t0, lambda, lam0;
2612 gmx_bool bSparse; /* use sparse matrix storage format */
2614 gmx_sparsematrix_t * sparse_matrix = NULL;
2615 real * full_matrix = NULL;
2616 em_state_t * state_work;
2618 /* added with respect to mdrun */
2619 int i, j, k, row, col;
2620 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2626 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2629 state_work = init_em_state();
2631 /* Init em and store the local state in state_minimum */
2632 init_em(fplog, NM, cr, inputrec,
2633 state_global, top_global, state_work, &top,
2635 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2636 nfile, fnm, &outf, NULL);
2638 natoms = top_global->natoms;
2646 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2647 " which MIGHT not be accurate enough for normal mode analysis.\n"
2648 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2649 " are fairly modest even if you recompile in double precision.\n\n");
2653 /* Check if we can/should use sparse storage format.
2655 * Sparse format is only useful when the Hessian itself is sparse, which it
2656 * will be when we use a cutoff.
2657 * For small systems (n<1000) it is easier to always use full matrix format, though.
2659 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2661 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2664 else if (top_global->natoms < 1000)
2666 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2671 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2677 sz = DIM*top_global->natoms;
2679 fprintf(stderr, "Allocating Hessian memory...\n\n");
2683 sparse_matrix = gmx_sparsematrix_init(sz);
2684 sparse_matrix->compressed_symmetric = TRUE;
2688 snew(full_matrix, sz*sz);
2692 /* Initial values */
2693 t0 = inputrec->init_t;
2694 lam0 = inputrec->fepvals->init_lambda;
2702 /* Write start time and temperature */
2703 print_em_start(fplog, cr, runtime, wcycle, NM);
2705 /* fudge nr of steps to nr of atoms */
2706 inputrec->nsteps = natoms*2;
2710 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2711 *(top_global->name), (int)inputrec->nsteps);
2714 nnodes = cr->nnodes;
2716 /* Make evaluate_energy do a single node force calculation */
2718 evaluate_energy(fplog, bVerbose, cr,
2719 state_global, top_global, state_work, top,
2720 inputrec, nrnb, wcycle, gstat,
2721 vsite, constr, fcd, graph, mdatoms, fr,
2722 mu_tot, enerd, vir, pres, -1, TRUE);
2723 cr->nnodes = nnodes;
2725 /* if forces are not small, warn user */
2726 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2728 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2729 if (state_work->fmax > 1.0e-3)
2731 md_print_info(cr, fplog,
2732 "The force is probably not small enough to "
2733 "ensure that you are at a minimum.\n"
2734 "Be aware that negative eigenvalues may occur\n"
2735 "when the resulting matrix is diagonalized.\n\n");
2738 /***********************************************************
2740 * Loop over all pairs in matrix
2742 * do_force called twice. Once with positive and
2743 * once with negative displacement
2745 ************************************************************/
2747 /* Steps are divided one by one over the nodes */
2748 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2751 for (d = 0; d < DIM; d++)
2753 x_min = state_work->s.x[atom][d];
2755 state_work->s.x[atom][d] = x_min - der_range;
2757 /* Make evaluate_energy do a single node force calculation */
2759 evaluate_energy(fplog, bVerbose, cr,
2760 state_global, top_global, state_work, top,
2761 inputrec, nrnb, wcycle, gstat,
2762 vsite, constr, fcd, graph, mdatoms, fr,
2763 mu_tot, enerd, vir, pres, atom*2, FALSE);
2765 for (i = 0; i < natoms; i++)
2767 copy_rvec(state_work->f[i], fneg[i]);
2770 state_work->s.x[atom][d] = x_min + der_range;
2772 evaluate_energy(fplog, bVerbose, cr,
2773 state_global, top_global, state_work, top,
2774 inputrec, nrnb, wcycle, gstat,
2775 vsite, constr, fcd, graph, mdatoms, fr,
2776 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2777 cr->nnodes = nnodes;
2779 /* x is restored to original */
2780 state_work->s.x[atom][d] = x_min;
2782 for (j = 0; j < natoms; j++)
2784 for (k = 0; (k < DIM); k++)
2787 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2795 #define mpi_type MPI_DOUBLE
2797 #define mpi_type MPI_FLOAT
2799 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2800 cr->mpi_comm_mygroup);
2805 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2811 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2812 cr->mpi_comm_mygroup, &stat);
2817 row = (atom + node)*DIM + d;
2819 for (j = 0; j < natoms; j++)
2821 for (k = 0; k < DIM; k++)
2827 if (col >= row && dfdx[j][k] != 0.0)
2829 gmx_sparsematrix_increment_value(sparse_matrix,
2830 row, col, dfdx[j][k]);
2835 full_matrix[row*sz+col] = dfdx[j][k];
2842 if (bVerbose && fplog)
2847 /* write progress */
2848 if (MASTER(cr) && bVerbose)
2850 fprintf(stderr, "\rFinished step %d out of %d",
2851 min(atom+nnodes, natoms), natoms);
2858 fprintf(stderr, "\n\nWriting Hessian...\n");
2859 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2862 finish_em(fplog, cr, outf, runtime, wcycle);
2864 runtime->nsteps_done = natoms*2;