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 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
52 #include "gromacs/random/random.h"
54 #include "gmx_fatal.h"
65 #include "md_support.h"
70 #include "mtop_util.h"
73 #include "gmx_omp_nthreads.h"
74 #include "md_logging.h"
76 #include "gromacs/fileio/confio.h"
77 #include "gromacs/fileio/trajectory_writing.h"
78 #include "gromacs/linearalgebra/mtxio.h"
79 #include "gromacs/linearalgebra/sparsematrix.h"
80 #include "gromacs/timing/wallcycle.h"
81 #include "gromacs/timing/walltime_accounting.h"
92 static em_state_t *init_em_state()
98 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
99 snew(ems->s.lambda, efptNR);
104 static void print_em_start(FILE *fplog,
106 gmx_walltime_accounting_t walltime_accounting,
107 gmx_wallcycle_t wcycle,
110 walltime_accounting_start(walltime_accounting);
111 wallcycle_start(wcycle, ewcRUN);
112 print_start(fplog, cr, walltime_accounting, name);
114 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
115 gmx_wallcycle_t wcycle)
117 wallcycle_stop(wcycle, ewcRUN);
119 walltime_accounting_end(walltime_accounting);
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)
136 "\nEnergy minimization reached the maximum number "
137 "of steps before the forces reached the requested "
138 "precision Fmax < %g.\n", ftol);
143 "\nEnergy minimization has stopped, but the forces have "
144 "not converged to the requested precision Fmax < %g (which "
145 "may not be possible for your system). It stopped "
146 "because the algorithm tried to make a new step whose size "
147 "was too small, or there was no change in the energy since "
148 "last step. Either way, we regard the minimization as "
149 "converged to within the available machine precision, "
150 "given your starting configuration and EM parameters.\n%s%s",
152 sizeof(real) < sizeof(double) ?
153 "\nDouble precision normally gives you higher accuracy, but "
154 "this is often not needed for preparing to run molecular "
158 "You might need to increase your constraint accuracy, or turn\n"
159 "off constraints altogether (set constraints = none in mdp file)\n" :
162 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
167 static void print_converged(FILE *fp, const char *alg, real ftol,
168 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
169 real epot, real fmax, int nfmax, real fnorm)
171 char buf[STEPSTRSIZE];
175 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
176 alg, ftol, gmx_step_str(count, buf));
178 else if (count < nsteps)
180 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
181 "but did not reach the requested Fmax < %g.\n",
182 alg, gmx_step_str(count, buf), ftol);
186 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
187 alg, ftol, gmx_step_str(count, buf));
191 fprintf(fp, "Potential Energy = %21.14e\n", epot);
192 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
193 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
195 fprintf(fp, "Potential Energy = %14.7e\n", epot);
196 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
197 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
201 static void get_f_norm_max(t_commrec *cr,
202 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
203 real *fnorm, real *fmax, int *a_fmax)
206 real fmax2, fmax2_0, fam;
207 int la_max, a_max, start, end, i, m, gf;
209 /* This routine finds the largest force and returns it.
210 * On parallel machines the global max is taken.
217 end = mdatoms->homenr;
218 if (mdatoms->cFREEZE)
220 for (i = start; i < end; i++)
222 gf = mdatoms->cFREEZE[i];
224 for (m = 0; m < DIM; m++)
226 if (!opts->nFreeze[gf][m])
241 for (i = start; i < end; i++)
253 if (la_max >= 0 && DOMAINDECOMP(cr))
255 a_max = cr->dd->gatindex[la_max];
263 snew(sum, 2*cr->nnodes+1);
264 sum[2*cr->nodeid] = fmax2;
265 sum[2*cr->nodeid+1] = a_max;
266 sum[2*cr->nnodes] = fnorm2;
267 gmx_sumd(2*cr->nnodes+1, sum, cr);
268 fnorm2 = sum[2*cr->nnodes];
269 /* Determine the global maximum */
270 for (i = 0; i < cr->nnodes; i++)
272 if (sum[2*i] > fmax2)
275 a_max = (int)(sum[2*i+1] + 0.5);
283 *fnorm = sqrt(fnorm2);
295 static void get_state_f_norm_max(t_commrec *cr,
296 t_grpopts *opts, t_mdatoms *mdatoms,
299 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
302 void init_em(FILE *fplog, const char *title,
303 t_commrec *cr, t_inputrec *ir,
304 t_state *state_global, gmx_mtop_t *top_global,
305 em_state_t *ems, gmx_localtop_t **top,
306 rvec **f, rvec **f_global,
307 t_nrnb *nrnb, rvec mu_tot,
308 t_forcerec *fr, gmx_enerdata_t **enerd,
309 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
310 gmx_vsite_t *vsite, gmx_constr_t constr,
311 int nfile, const t_filenm fnm[],
312 gmx_mdoutf_t *outf, t_mdebin **mdebin)
319 fprintf(fplog, "Initiating %s\n", title);
322 state_global->ngtc = 0;
324 /* Initialize lambda variables */
325 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
329 if (DOMAINDECOMP(cr))
331 *top = dd_init_local_top(top_global);
333 dd_init_local_state(cr->dd, state_global, &ems->s);
337 /* Distribute the charge groups over the nodes from the master node */
338 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
339 state_global, top_global, ir,
340 &ems->s, &ems->f, mdatoms, *top,
341 fr, vsite, NULL, constr,
343 dd_store_state(cr->dd, &ems->s);
347 snew(*f_global, top_global->natoms);
357 snew(*f, top_global->natoms);
359 /* Just copy the state */
360 ems->s = *state_global;
361 snew(ems->s.x, ems->s.nalloc);
362 snew(ems->f, ems->s.nalloc);
363 for (i = 0; i < state_global->natoms; i++)
365 copy_rvec(state_global->x[i], ems->s.x[i]);
367 copy_mat(state_global->box, ems->s.box);
369 *top = gmx_mtop_generate_local_top(top_global, ir);
372 forcerec_set_excl_load(fr, *top);
374 setup_bonded_threading(fr, &(*top)->idef);
376 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
378 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
385 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
386 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
390 set_vsite_top(vsite, *top, mdatoms, cr);
396 if (ir->eConstrAlg == econtSHAKE &&
397 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
399 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
400 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
403 if (!DOMAINDECOMP(cr))
405 set_constraints(constr, *top, ir, mdatoms, cr);
408 if (!ir->bContinuation)
410 /* Constrain the starting coordinates */
412 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
413 ir, NULL, cr, -1, 0, mdatoms,
414 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
415 ems->s.lambda[efptFEP], &dvdl_constr,
416 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
422 *gstat = global_stat_init(ir);
425 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, NULL);
428 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
433 /* Init bin for energy stuff */
434 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
438 calc_shifts(ems->s.box, fr->shift_vec);
441 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
442 gmx_walltime_accounting_t walltime_accounting,
443 gmx_wallcycle_t wcycle)
445 if (!(cr->duty & DUTY_PME))
447 /* Tell the PME only node to finish */
448 gmx_pme_send_finish(cr);
453 em_time_end(walltime_accounting, wcycle);
456 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
465 static void copy_em_coords(em_state_t *ems, t_state *state)
469 for (i = 0; (i < state->natoms); i++)
471 copy_rvec(ems->s.x[i], state->x[i]);
475 static void write_em_traj(FILE *fplog, t_commrec *cr,
477 gmx_bool bX, gmx_bool bF, const char *confout,
478 gmx_mtop_t *top_global,
479 t_inputrec *ir, gmx_int64_t step,
481 t_state *state_global, rvec *f_global)
485 if ((bX || bF || confout != NULL) && !DOMAINDECOMP(cr))
487 copy_em_coords(state, state_global);
494 mdof_flags |= MDOF_X;
498 mdof_flags |= MDOF_F;
500 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
501 top_global, step, (double)step,
502 &state->s, state_global, state->f, f_global);
504 if (confout != NULL && MASTER(cr))
506 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
508 /* Make molecules whole only for confout writing */
509 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
513 write_sto_conf_mtop(confout,
514 *top_global->name, top_global,
515 state_global->x, NULL, ir->ePBC, state_global->box);
519 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
521 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
522 gmx_constr_t constr, gmx_localtop_t *top,
523 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
536 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
538 gmx_incons("state mismatch in do_em_step");
541 s2->flags = s1->flags;
543 if (s2->nalloc != s1->nalloc)
545 s2->nalloc = s1->nalloc;
546 srenew(s2->x, s1->nalloc);
547 srenew(ems2->f, s1->nalloc);
548 if (s2->flags & (1<<estCGP))
550 srenew(s2->cg_p, s1->nalloc);
554 s2->natoms = s1->natoms;
555 copy_mat(s1->box, s2->box);
556 /* Copy free energy state */
557 for (i = 0; i < efptNR; i++)
559 s2->lambda[i] = s1->lambda[i];
561 copy_mat(s1->box, s2->box);
569 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
574 #pragma omp for schedule(static) nowait
575 for (i = start; i < end; i++)
581 for (m = 0; m < DIM; m++)
583 if (ir->opts.nFreeze[gf][m])
589 x2[i][m] = x1[i][m] + a*f[i][m];
594 if (s2->flags & (1<<estCGP))
596 /* Copy the CG p vector */
599 #pragma omp for schedule(static) nowait
600 for (i = start; i < end; i++)
602 copy_rvec(x1[i], x2[i]);
606 if (DOMAINDECOMP(cr))
608 s2->ddp_count = s1->ddp_count;
609 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
612 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
613 srenew(s2->cg_gl, s2->cg_gl_nalloc);
616 s2->ncg_gl = s1->ncg_gl;
617 #pragma omp for schedule(static) nowait
618 for (i = 0; i < s2->ncg_gl; i++)
620 s2->cg_gl[i] = s1->cg_gl[i];
622 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
628 wallcycle_start(wcycle, ewcCONSTR);
630 constrain(NULL, TRUE, TRUE, constr, &top->idef,
631 ir, NULL, cr, count, 0, md,
632 s1->x, s2->x, NULL, bMolPBC, s2->box,
633 s2->lambda[efptBONDED], &dvdl_constr,
634 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
635 wallcycle_stop(wcycle, ewcCONSTR);
639 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
640 gmx_mtop_t *top_global, t_inputrec *ir,
641 em_state_t *ems, gmx_localtop_t *top,
642 t_mdatoms *mdatoms, t_forcerec *fr,
643 gmx_vsite_t *vsite, gmx_constr_t constr,
644 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
646 /* Repartition the domain decomposition */
647 wallcycle_start(wcycle, ewcDOMDEC);
648 dd_partition_system(fplog, step, cr, FALSE, 1,
649 NULL, top_global, ir,
651 mdatoms, top, fr, vsite, NULL, constr,
652 nrnb, wcycle, FALSE);
653 dd_store_state(cr->dd, &ems->s);
654 wallcycle_stop(wcycle, ewcDOMDEC);
657 static void evaluate_energy(FILE *fplog, t_commrec *cr,
658 gmx_mtop_t *top_global,
659 em_state_t *ems, gmx_localtop_t *top,
660 t_inputrec *inputrec,
661 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
662 gmx_global_stat_t gstat,
663 gmx_vsite_t *vsite, gmx_constr_t constr,
665 t_graph *graph, t_mdatoms *mdatoms,
666 t_forcerec *fr, rvec mu_tot,
667 gmx_enerdata_t *enerd, tensor vir, tensor pres,
668 gmx_int64_t count, gmx_bool bFirst)
673 tensor force_vir, shake_vir, ekin;
674 real dvdl_constr, prescorr, enercorr, dvdlcorr;
677 /* Set the time to the initial time, the time does not change during EM */
678 t = inputrec->init_t;
681 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
683 /* This is the first state or an old state used before the last ns */
689 if (inputrec->nstlist > 0)
693 else if (inputrec->nstlist == -1)
695 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
698 gmx_sumi(1, &nabnsb, cr);
706 construct_vsites(vsite, ems->s.x, 1, NULL,
707 top->idef.iparams, top->idef.il,
708 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
711 if (DOMAINDECOMP(cr) && bNS)
713 /* Repartition the domain decomposition */
714 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
715 ems, top, mdatoms, fr, vsite, constr,
719 /* Calc force & energy on new trial position */
720 /* do_force always puts the charge groups in the box and shifts again
721 * We do not unshift, so molecules are always whole in congrad.c
723 do_force(fplog, cr, inputrec,
724 count, nrnb, wcycle, top, &top_global->groups,
725 ems->s.box, ems->s.x, &ems->s.hist,
726 ems->f, force_vir, mdatoms, enerd, fcd,
727 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
728 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
729 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
730 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
732 /* Clear the unused shake virial and pressure */
733 clear_mat(shake_vir);
736 /* Communicate stuff when parallel */
737 if (PAR(cr) && inputrec->eI != eiNM)
739 wallcycle_start(wcycle, ewcMoveE);
741 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
742 inputrec, NULL, NULL, NULL, 1, &terminate,
743 top_global, &ems->s, FALSE,
749 wallcycle_stop(wcycle, ewcMoveE);
752 /* Calculate long range corrections to pressure and energy */
753 calc_dispcorr(fplog, inputrec, fr, count, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
754 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
755 enerd->term[F_DISPCORR] = enercorr;
756 enerd->term[F_EPOT] += enercorr;
757 enerd->term[F_PRES] += prescorr;
758 enerd->term[F_DVDL] += dvdlcorr;
760 ems->epot = enerd->term[F_EPOT];
764 /* Project out the constraint components of the force */
765 wallcycle_start(wcycle, ewcCONSTR);
767 constrain(NULL, FALSE, FALSE, constr, &top->idef,
768 inputrec, NULL, cr, count, 0, mdatoms,
769 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
770 ems->s.lambda[efptBONDED], &dvdl_constr,
771 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
772 if (fr->bSepDVDL && fplog)
774 gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
776 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
777 m_add(force_vir, shake_vir, vir);
778 wallcycle_stop(wcycle, ewcCONSTR);
782 copy_mat(force_vir, vir);
786 enerd->term[F_PRES] =
787 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
789 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
791 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
793 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
797 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
799 em_state_t *s_min, em_state_t *s_b)
803 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
805 unsigned char *grpnrFREEZE;
809 fprintf(debug, "Doing reorder_partsum\n");
815 cgs_gl = dd_charge_groups_global(cr->dd);
816 index = cgs_gl->index;
818 /* Collect fm in a global vector fmg.
819 * This conflicts with the spirit of domain decomposition,
820 * but to fully optimize this a much more complicated algorithm is required.
822 snew(fmg, mtop->natoms);
824 ncg = s_min->s.ncg_gl;
825 cg_gl = s_min->s.cg_gl;
827 for (c = 0; c < ncg; c++)
832 for (a = a0; a < a1; a++)
834 copy_rvec(fm[i], fmg[a]);
838 gmx_sum(mtop->natoms*3, fmg[0], cr);
840 /* Now we will determine the part of the sum for the cgs in state s_b */
842 cg_gl = s_b->s.cg_gl;
846 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
847 for (c = 0; c < ncg; c++)
852 for (a = a0; a < a1; a++)
854 if (mdatoms->cFREEZE && grpnrFREEZE)
858 for (m = 0; m < DIM; m++)
860 if (!opts->nFreeze[gf][m])
862 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
874 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
876 em_state_t *s_min, em_state_t *s_b)
882 /* This is just the classical Polak-Ribiere calculation of beta;
883 * it looks a bit complicated since we take freeze groups into account,
884 * and might have to sum it in parallel runs.
887 if (!DOMAINDECOMP(cr) ||
888 (s_min->s.ddp_count == cr->dd->ddp_count &&
889 s_b->s.ddp_count == cr->dd->ddp_count))
895 /* This part of code can be incorrect with DD,
896 * since the atom ordering in s_b and s_min might differ.
898 for (i = 0; i < mdatoms->homenr; i++)
900 if (mdatoms->cFREEZE)
902 gf = mdatoms->cFREEZE[i];
904 for (m = 0; m < DIM; m++)
906 if (!opts->nFreeze[gf][m])
908 sum += (fb[i][m] - fm[i][m])*fb[i][m];
915 /* We need to reorder cgs while summing */
916 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
920 gmx_sumd(1, &sum, cr);
923 return sum/sqr(s_min->fnorm);
926 double do_cg(FILE *fplog, t_commrec *cr,
927 int nfile, const t_filenm fnm[],
928 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
929 int gmx_unused nstglobalcomm,
930 gmx_vsite_t *vsite, gmx_constr_t constr,
931 int gmx_unused stepout,
932 t_inputrec *inputrec,
933 gmx_mtop_t *top_global, t_fcdata *fcd,
934 t_state *state_global,
936 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
937 gmx_edsam_t gmx_unused ed,
939 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
940 gmx_membed_t gmx_unused membed,
941 real gmx_unused cpt_period, real gmx_unused max_hours,
942 const char gmx_unused *deviceOptions,
943 unsigned long gmx_unused Flags,
944 gmx_walltime_accounting_t walltime_accounting)
946 const char *CG = "Polak-Ribiere Conjugate Gradients";
948 em_state_t *s_min, *s_a, *s_b, *s_c;
950 gmx_enerdata_t *enerd;
952 gmx_global_stat_t gstat;
954 rvec *f_global, *p, *sf, *sfm;
955 double gpa, gpb, gpc, tmp, sum[2], minstep;
958 real a, b, c, beta = 0.0;
962 gmx_bool converged, foundlower;
964 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
966 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
968 int i, m, gf, step, nminstep;
973 s_min = init_em_state();
974 s_a = init_em_state();
975 s_b = init_em_state();
976 s_c = init_em_state();
978 /* Init em and store the local state in s_min */
979 init_em(fplog, CG, cr, inputrec,
980 state_global, top_global, s_min, &top, &f, &f_global,
981 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
982 nfile, fnm, &outf, &mdebin);
984 /* Print to log file */
985 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
987 /* Max number of steps */
988 number_steps = inputrec->nsteps;
992 sp_header(stderr, CG, inputrec->em_tol, number_steps);
996 sp_header(fplog, CG, inputrec->em_tol, number_steps);
999 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1000 /* do_force always puts the charge groups in the box and shifts again
1001 * We do not unshift, so molecules are always whole in congrad.c
1003 evaluate_energy(fplog, cr,
1004 top_global, s_min, top,
1005 inputrec, nrnb, wcycle, gstat,
1006 vsite, constr, fcd, graph, mdatoms, fr,
1007 mu_tot, enerd, vir, pres, -1, TRUE);
1012 /* Copy stuff to the energy bin for easy printing etc. */
1013 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1014 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1015 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1017 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1018 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1019 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1023 /* Estimate/guess the initial stepsize */
1024 stepsize = inputrec->em_stepsize/s_min->fnorm;
1028 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1029 s_min->fmax, s_min->a_fmax+1);
1030 fprintf(stderr, " F-Norm = %12.5e\n",
1031 s_min->fnorm/sqrt(state_global->natoms));
1032 fprintf(stderr, "\n");
1033 /* and copy to the log file too... */
1034 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1035 s_min->fmax, s_min->a_fmax+1);
1036 fprintf(fplog, " F-Norm = %12.5e\n",
1037 s_min->fnorm/sqrt(state_global->natoms));
1038 fprintf(fplog, "\n");
1040 /* Start the loop over CG steps.
1041 * Each successful step is counted, and we continue until
1042 * we either converge or reach the max number of steps.
1045 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1048 /* start taking steps in a new direction
1049 * First time we enter the routine, beta=0, and the direction is
1050 * simply the negative gradient.
1053 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1058 for (i = 0; i < mdatoms->homenr; i++)
1060 if (mdatoms->cFREEZE)
1062 gf = mdatoms->cFREEZE[i];
1064 for (m = 0; m < DIM; m++)
1066 if (!inputrec->opts.nFreeze[gf][m])
1068 p[i][m] = sf[i][m] + beta*p[i][m];
1069 gpa -= p[i][m]*sf[i][m];
1070 /* f is negative gradient, thus the sign */
1079 /* Sum the gradient along the line across CPUs */
1082 gmx_sumd(1, &gpa, cr);
1085 /* Calculate the norm of the search vector */
1086 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1088 /* Just in case stepsize reaches zero due to numerical precision... */
1091 stepsize = inputrec->em_stepsize/pnorm;
1095 * Double check the value of the derivative in the search direction.
1096 * If it is positive it must be due to the old information in the
1097 * CG formula, so just remove that and start over with beta=0.
1098 * This corresponds to a steepest descent step.
1103 step--; /* Don't count this step since we are restarting */
1104 continue; /* Go back to the beginning of the big for-loop */
1107 /* Calculate minimum allowed stepsize, before the average (norm)
1108 * relative change in coordinate is smaller than precision
1111 for (i = 0; i < mdatoms->homenr; i++)
1113 for (m = 0; m < DIM; m++)
1115 tmp = fabs(s_min->s.x[i][m]);
1124 /* Add up from all CPUs */
1127 gmx_sumd(1, &minstep, cr);
1130 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1132 if (stepsize < minstep)
1138 /* Write coordinates if necessary */
1139 do_x = do_per_step(step, inputrec->nstxout);
1140 do_f = do_per_step(step, inputrec->nstfout);
1142 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1143 top_global, inputrec, step,
1144 s_min, state_global, f_global);
1146 /* Take a step downhill.
1147 * In theory, we should minimize the function along this direction.
1148 * That is quite possible, but it turns out to take 5-10 function evaluations
1149 * for each line. However, we dont really need to find the exact minimum -
1150 * it is much better to start a new CG step in a modified direction as soon
1151 * as we are close to it. This will save a lot of energy evaluations.
1153 * In practice, we just try to take a single step.
1154 * If it worked (i.e. lowered the energy), we increase the stepsize but
1155 * the continue straight to the next CG step without trying to find any minimum.
1156 * If it didn't work (higher energy), there must be a minimum somewhere between
1157 * the old position and the new one.
1159 * Due to the finite numerical accuracy, it turns out that it is a good idea
1160 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1161 * This leads to lower final energies in the tests I've done. / Erik
1163 s_a->epot = s_min->epot;
1165 c = a + stepsize; /* reference position along line is zero */
1167 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1169 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1170 s_min, top, mdatoms, fr, vsite, constr,
1174 /* Take a trial step (new coords in s_c) */
1175 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1176 constr, top, nrnb, wcycle, -1);
1179 /* Calculate energy for the trial step */
1180 evaluate_energy(fplog, cr,
1181 top_global, s_c, top,
1182 inputrec, nrnb, wcycle, gstat,
1183 vsite, constr, fcd, graph, mdatoms, fr,
1184 mu_tot, enerd, vir, pres, -1, FALSE);
1186 /* Calc derivative along line */
1190 for (i = 0; i < mdatoms->homenr; i++)
1192 for (m = 0; m < DIM; m++)
1194 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1197 /* Sum the gradient along the line across CPUs */
1200 gmx_sumd(1, &gpc, cr);
1203 /* This is the max amount of increase in energy we tolerate */
1204 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1206 /* Accept the step if the energy is lower, or if it is not significantly higher
1207 * and the line derivative is still negative.
1209 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1212 /* Great, we found a better energy. Increase step for next iteration
1213 * if we are still going down, decrease it otherwise
1217 stepsize *= 1.618034; /* The golden section */
1221 stepsize *= 0.618034; /* 1/golden section */
1226 /* New energy is the same or higher. We will have to do some work
1227 * to find a smaller value in the interval. Take smaller step next time!
1230 stepsize *= 0.618034;
1236 /* OK, if we didn't find a lower value we will have to locate one now - there must
1237 * be one in the interval [a=0,c].
1238 * The same thing is valid here, though: Don't spend dozens of iterations to find
1239 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1240 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1242 * I also have a safeguard for potentially really patological functions so we never
1243 * take more than 20 steps before we give up ...
1245 * If we already found a lower value we just skip this step and continue to the update.
1253 /* Select a new trial point.
1254 * If the derivatives at points a & c have different sign we interpolate to zero,
1255 * otherwise just do a bisection.
1257 if (gpa < 0 && gpc > 0)
1259 b = a + gpa*(a-c)/(gpc-gpa);
1266 /* safeguard if interpolation close to machine accuracy causes errors:
1267 * never go outside the interval
1269 if (b <= a || b >= c)
1274 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1276 /* Reload the old state */
1277 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1278 s_min, top, mdatoms, fr, vsite, constr,
1282 /* Take a trial step to this new point - new coords in s_b */
1283 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1284 constr, top, nrnb, wcycle, -1);
1287 /* Calculate energy for the trial step */
1288 evaluate_energy(fplog, cr,
1289 top_global, s_b, top,
1290 inputrec, nrnb, wcycle, gstat,
1291 vsite, constr, fcd, graph, mdatoms, fr,
1292 mu_tot, enerd, vir, pres, -1, FALSE);
1294 /* p does not change within a step, but since the domain decomposition
1295 * might change, we have to use cg_p of s_b here.
1300 for (i = 0; i < mdatoms->homenr; i++)
1302 for (m = 0; m < DIM; m++)
1304 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1307 /* Sum the gradient along the line across CPUs */
1310 gmx_sumd(1, &gpb, cr);
1315 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1316 s_a->epot, s_b->epot, s_c->epot, gpb);
1319 epot_repl = s_b->epot;
1321 /* Keep one of the intervals based on the value of the derivative at the new point */
1324 /* Replace c endpoint with b */
1325 swap_em_state(s_b, s_c);
1331 /* Replace a endpoint with b */
1332 swap_em_state(s_b, s_a);
1338 * Stop search as soon as we find a value smaller than the endpoints.
1339 * Never run more than 20 steps, no matter what.
1343 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1346 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1349 /* OK. We couldn't find a significantly lower energy.
1350 * If beta==0 this was steepest descent, and then we give up.
1351 * If not, set beta=0 and restart with steepest descent before quitting.
1361 /* Reset memory before giving up */
1367 /* Select min energy state of A & C, put the best in B.
1369 if (s_c->epot < s_a->epot)
1373 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1374 s_c->epot, s_a->epot);
1376 swap_em_state(s_b, s_c);
1384 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1385 s_a->epot, s_c->epot);
1387 swap_em_state(s_b, s_a);
1397 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1400 swap_em_state(s_b, s_c);
1405 /* new search direction */
1406 /* beta = 0 means forget all memory and restart with steepest descents. */
1407 if (nstcg && ((step % nstcg) == 0))
1413 /* s_min->fnorm cannot be zero, because then we would have converged
1417 /* Polak-Ribiere update.
1418 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1420 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1422 /* Limit beta to prevent oscillations */
1423 if (fabs(beta) > 5.0)
1429 /* update positions */
1430 swap_em_state(s_min, s_b);
1433 /* Print it if necessary */
1438 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1439 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1440 s_min->fmax, s_min->a_fmax+1);
1442 /* Store the new (lower) energies */
1443 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1444 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1445 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1447 do_log = do_per_step(step, inputrec->nstlog);
1448 do_ene = do_per_step(step, inputrec->nstenergy);
1451 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1453 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1454 do_log ? fplog : NULL, step, step, eprNORMAL,
1455 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1458 /* Stop when the maximum force lies below tolerance.
1459 * If we have reached machine precision, converged is already set to true.
1461 converged = converged || (s_min->fmax < inputrec->em_tol);
1463 } /* End of the loop */
1467 step--; /* we never took that last step in this case */
1470 if (s_min->fmax > inputrec->em_tol)
1474 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1475 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1482 /* If we printed energy and/or logfile last step (which was the last step)
1483 * we don't have to do it again, but otherwise print the final values.
1487 /* Write final value to log since we didn't do anything the last step */
1488 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1490 if (!do_ene || !do_log)
1492 /* Write final energy file entries */
1493 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1494 !do_log ? fplog : NULL, step, step, eprNORMAL,
1495 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1499 /* Print some stuff... */
1502 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1506 * For accurate normal mode calculation it is imperative that we
1507 * store the last conformation into the full precision binary trajectory.
1509 * However, we should only do it if we did NOT already write this step
1510 * above (which we did if do_x or do_f was true).
1512 do_x = !do_per_step(step, inputrec->nstxout);
1513 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1515 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1516 top_global, inputrec, step,
1517 s_min, state_global, f_global);
1519 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1523 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1524 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1525 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1526 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1528 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1531 finish_em(cr, outf, walltime_accounting, wcycle);
1533 /* To print the actual number of steps we needed somewhere */
1534 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1537 } /* That's all folks */
1540 double do_lbfgs(FILE *fplog, t_commrec *cr,
1541 int nfile, const t_filenm fnm[],
1542 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1543 int gmx_unused nstglobalcomm,
1544 gmx_vsite_t *vsite, gmx_constr_t constr,
1545 int gmx_unused stepout,
1546 t_inputrec *inputrec,
1547 gmx_mtop_t *top_global, t_fcdata *fcd,
1550 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1551 gmx_edsam_t gmx_unused ed,
1553 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1554 gmx_membed_t gmx_unused membed,
1555 real gmx_unused cpt_period, real gmx_unused max_hours,
1556 const char gmx_unused *deviceOptions,
1557 unsigned long gmx_unused Flags,
1558 gmx_walltime_accounting_t walltime_accounting)
1560 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1562 gmx_localtop_t *top;
1563 gmx_enerdata_t *enerd;
1565 gmx_global_stat_t gstat;
1568 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1569 double stepsize, gpa, gpb, gpc, tmp, minstep;
1570 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1571 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1572 real a, b, c, maxdelta, delta;
1573 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1574 real dgdx, dgdg, sq, yr, beta;
1576 gmx_bool converged, first;
1579 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1581 int start, end, number_steps;
1583 int i, k, m, n, nfmax, gf, step;
1590 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1595 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).");
1598 n = 3*state->natoms;
1599 nmaxcorr = inputrec->nbfgscorr;
1601 /* Allocate memory */
1602 /* Use pointers to real so we dont have to loop over both atoms and
1603 * dimensions all the time...
1604 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1605 * that point to the same memory.
1618 snew(rho, nmaxcorr);
1619 snew(alpha, nmaxcorr);
1622 for (i = 0; i < nmaxcorr; i++)
1628 for (i = 0; i < nmaxcorr; i++)
1637 init_em(fplog, LBFGS, cr, inputrec,
1638 state, top_global, &ems, &top, &f, &f_global,
1639 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1640 nfile, fnm, &outf, &mdebin);
1641 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1642 * so we free some memory again.
1647 xx = (real *)state->x;
1651 end = mdatoms->homenr;
1653 /* Print to log file */
1654 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1656 do_log = do_ene = do_x = do_f = TRUE;
1658 /* Max number of steps */
1659 number_steps = inputrec->nsteps;
1661 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1663 for (i = start; i < end; i++)
1665 if (mdatoms->cFREEZE)
1667 gf = mdatoms->cFREEZE[i];
1669 for (m = 0; m < DIM; m++)
1671 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1676 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1680 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1685 construct_vsites(vsite, state->x, 1, NULL,
1686 top->idef.iparams, top->idef.il,
1687 fr->ePBC, fr->bMolPBC, cr, state->box);
1690 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1691 /* do_force always puts the charge groups in the box and shifts again
1692 * We do not unshift, so molecules are always whole
1697 evaluate_energy(fplog, cr,
1698 top_global, &ems, top,
1699 inputrec, nrnb, wcycle, gstat,
1700 vsite, constr, fcd, graph, mdatoms, fr,
1701 mu_tot, enerd, vir, pres, -1, TRUE);
1706 /* Copy stuff to the energy bin for easy printing etc. */
1707 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1708 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1709 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1711 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1712 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1713 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1717 /* This is the starting energy */
1718 Epot = enerd->term[F_EPOT];
1724 /* Set the initial step.
1725 * since it will be multiplied by the non-normalized search direction
1726 * vector (force vector the first time), we scale it by the
1727 * norm of the force.
1732 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1733 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1734 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1735 fprintf(stderr, "\n");
1736 /* and copy to the log file too... */
1737 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1738 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1739 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1740 fprintf(fplog, "\n");
1744 for (i = 0; i < n; i++)
1748 dx[point][i] = ff[i]; /* Initial search direction */
1756 stepsize = 1.0/fnorm;
1759 /* Start the loop over BFGS steps.
1760 * Each successful step is counted, and we continue until
1761 * we either converge or reach the max number of steps.
1766 /* Set the gradient from the force */
1768 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1771 /* Write coordinates if necessary */
1772 do_x = do_per_step(step, inputrec->nstxout);
1773 do_f = do_per_step(step, inputrec->nstfout);
1778 mdof_flags |= MDOF_X;
1783 mdof_flags |= MDOF_F;
1786 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1787 top_global, step, (real)step, state, state, f, f);
1789 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1791 /* pointer to current direction - point=0 first time here */
1794 /* calculate line gradient */
1795 for (gpa = 0, i = 0; i < n; i++)
1800 /* Calculate minimum allowed stepsize, before the average (norm)
1801 * relative change in coordinate is smaller than precision
1803 for (minstep = 0, i = 0; i < n; i++)
1813 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1815 if (stepsize < minstep)
1821 /* Store old forces and coordinates */
1822 for (i = 0; i < n; i++)
1831 for (i = 0; i < n; i++)
1836 /* Take a step downhill.
1837 * In theory, we should minimize the function along this direction.
1838 * That is quite possible, but it turns out to take 5-10 function evaluations
1839 * for each line. However, we dont really need to find the exact minimum -
1840 * it is much better to start a new BFGS step in a modified direction as soon
1841 * as we are close to it. This will save a lot of energy evaluations.
1843 * In practice, we just try to take a single step.
1844 * If it worked (i.e. lowered the energy), we increase the stepsize but
1845 * the continue straight to the next BFGS step without trying to find any minimum.
1846 * If it didn't work (higher energy), there must be a minimum somewhere between
1847 * the old position and the new one.
1849 * Due to the finite numerical accuracy, it turns out that it is a good idea
1850 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1851 * This leads to lower final energies in the tests I've done. / Erik
1856 c = a + stepsize; /* reference position along line is zero */
1858 /* Check stepsize first. We do not allow displacements
1859 * larger than emstep.
1865 for (i = 0; i < n; i++)
1868 if (delta > maxdelta)
1873 if (maxdelta > inputrec->em_stepsize)
1878 while (maxdelta > inputrec->em_stepsize);
1880 /* Take a trial step */
1881 for (i = 0; i < n; i++)
1883 xc[i] = lastx[i] + c*s[i];
1887 /* Calculate energy for the trial step */
1888 ems.s.x = (rvec *)xc;
1890 evaluate_energy(fplog, cr,
1891 top_global, &ems, top,
1892 inputrec, nrnb, wcycle, gstat,
1893 vsite, constr, fcd, graph, mdatoms, fr,
1894 mu_tot, enerd, vir, pres, step, FALSE);
1897 /* Calc derivative along line */
1898 for (gpc = 0, i = 0; i < n; i++)
1900 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1902 /* Sum the gradient along the line across CPUs */
1905 gmx_sumd(1, &gpc, cr);
1908 /* This is the max amount of increase in energy we tolerate */
1909 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1911 /* Accept the step if the energy is lower, or if it is not significantly higher
1912 * and the line derivative is still negative.
1914 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1917 /* Great, we found a better energy. Increase step for next iteration
1918 * if we are still going down, decrease it otherwise
1922 stepsize *= 1.618034; /* The golden section */
1926 stepsize *= 0.618034; /* 1/golden section */
1931 /* New energy is the same or higher. We will have to do some work
1932 * to find a smaller value in the interval. Take smaller step next time!
1935 stepsize *= 0.618034;
1938 /* OK, if we didn't find a lower value we will have to locate one now - there must
1939 * be one in the interval [a=0,c].
1940 * The same thing is valid here, though: Don't spend dozens of iterations to find
1941 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1942 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1944 * I also have a safeguard for potentially really patological functions so we never
1945 * take more than 20 steps before we give up ...
1947 * If we already found a lower value we just skip this step and continue to the update.
1956 /* Select a new trial point.
1957 * If the derivatives at points a & c have different sign we interpolate to zero,
1958 * otherwise just do a bisection.
1961 if (gpa < 0 && gpc > 0)
1963 b = a + gpa*(a-c)/(gpc-gpa);
1970 /* safeguard if interpolation close to machine accuracy causes errors:
1971 * never go outside the interval
1973 if (b <= a || b >= c)
1978 /* Take a trial step */
1979 for (i = 0; i < n; i++)
1981 xb[i] = lastx[i] + b*s[i];
1985 /* Calculate energy for the trial step */
1986 ems.s.x = (rvec *)xb;
1988 evaluate_energy(fplog, cr,
1989 top_global, &ems, top,
1990 inputrec, nrnb, wcycle, gstat,
1991 vsite, constr, fcd, graph, mdatoms, fr,
1992 mu_tot, enerd, vir, pres, step, FALSE);
1997 for (gpb = 0, i = 0; i < n; i++)
1999 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2002 /* Sum the gradient along the line across CPUs */
2005 gmx_sumd(1, &gpb, cr);
2008 /* Keep one of the intervals based on the value of the derivative at the new point */
2011 /* Replace c endpoint with b */
2015 /* swap coord pointers b/c */
2025 /* Replace a endpoint with b */
2029 /* swap coord pointers a/b */
2039 * Stop search as soon as we find a value smaller than the endpoints,
2040 * or if the tolerance is below machine precision.
2041 * Never run more than 20 steps, no matter what.
2045 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2047 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2049 /* OK. We couldn't find a significantly lower energy.
2050 * If ncorr==0 this was steepest descent, and then we give up.
2051 * If not, reset memory to restart as steepest descent before quitting.
2063 /* Search in gradient direction */
2064 for (i = 0; i < n; i++)
2066 dx[point][i] = ff[i];
2068 /* Reset stepsize */
2069 stepsize = 1.0/fnorm;
2074 /* Select min energy state of A & C, put the best in xx/ff/Epot
2080 for (i = 0; i < n; i++)
2091 for (i = 0; i < n; i++)
2105 for (i = 0; i < n; i++)
2113 /* Update the memory information, and calculate a new
2114 * approximation of the inverse hessian
2117 /* Have new data in Epot, xx, ff */
2118 if (ncorr < nmaxcorr)
2123 for (i = 0; i < n; i++)
2125 dg[point][i] = lastf[i]-ff[i];
2126 dx[point][i] *= stepsize;
2131 for (i = 0; i < n; i++)
2133 dgdg += dg[point][i]*dg[point][i];
2134 dgdx += dg[point][i]*dx[point][i];
2139 rho[point] = 1.0/dgdx;
2142 if (point >= nmaxcorr)
2148 for (i = 0; i < n; i++)
2155 /* Recursive update. First go back over the memory points */
2156 for (k = 0; k < ncorr; k++)
2165 for (i = 0; i < n; i++)
2167 sq += dx[cp][i]*p[i];
2170 alpha[cp] = rho[cp]*sq;
2172 for (i = 0; i < n; i++)
2174 p[i] -= alpha[cp]*dg[cp][i];
2178 for (i = 0; i < n; i++)
2183 /* And then go forward again */
2184 for (k = 0; k < ncorr; k++)
2187 for (i = 0; i < n; i++)
2189 yr += p[i]*dg[cp][i];
2193 beta = alpha[cp]-beta;
2195 for (i = 0; i < n; i++)
2197 p[i] += beta*dx[cp][i];
2207 for (i = 0; i < n; i++)
2211 dx[point][i] = p[i];
2221 /* Test whether the convergence criterion is met */
2222 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2224 /* Print it if necessary */
2229 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2230 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2232 /* Store the new (lower) energies */
2233 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2234 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2235 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2236 do_log = do_per_step(step, inputrec->nstlog);
2237 do_ene = do_per_step(step, inputrec->nstenergy);
2240 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2242 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2243 do_log ? fplog : NULL, step, step, eprNORMAL,
2244 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2247 /* Stop when the maximum force lies below tolerance.
2248 * If we have reached machine precision, converged is already set to true.
2251 converged = converged || (fmax < inputrec->em_tol);
2253 } /* End of the loop */
2257 step--; /* we never took that last step in this case */
2260 if (fmax > inputrec->em_tol)
2264 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2265 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2270 /* If we printed energy and/or logfile last step (which was the last step)
2271 * we don't have to do it again, but otherwise print the final values.
2273 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2275 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2277 if (!do_ene || !do_log) /* Write final energy file entries */
2279 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2280 !do_log ? fplog : NULL, step, step, eprNORMAL,
2281 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2284 /* Print some stuff... */
2287 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2291 * For accurate normal mode calculation it is imperative that we
2292 * store the last conformation into the full precision binary trajectory.
2294 * However, we should only do it if we did NOT already write this step
2295 * above (which we did if do_x or do_f was true).
2297 do_x = !do_per_step(step, inputrec->nstxout);
2298 do_f = !do_per_step(step, inputrec->nstfout);
2299 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2300 top_global, inputrec, step,
2305 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2306 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2307 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2308 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2310 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2313 finish_em(cr, outf, walltime_accounting, wcycle);
2315 /* To print the actual number of steps we needed somewhere */
2316 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2319 } /* That's all folks */
2322 double do_steep(FILE *fplog, t_commrec *cr,
2323 int nfile, const t_filenm fnm[],
2324 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2325 int gmx_unused nstglobalcomm,
2326 gmx_vsite_t *vsite, gmx_constr_t constr,
2327 int gmx_unused stepout,
2328 t_inputrec *inputrec,
2329 gmx_mtop_t *top_global, t_fcdata *fcd,
2330 t_state *state_global,
2332 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2333 gmx_edsam_t gmx_unused ed,
2335 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2336 gmx_membed_t gmx_unused membed,
2337 real gmx_unused cpt_period, real gmx_unused max_hours,
2338 const char gmx_unused *deviceOptions,
2339 unsigned long gmx_unused Flags,
2340 gmx_walltime_accounting_t walltime_accounting)
2342 const char *SD = "Steepest Descents";
2343 em_state_t *s_min, *s_try;
2345 gmx_localtop_t *top;
2346 gmx_enerdata_t *enerd;
2348 gmx_global_stat_t gstat;
2350 real stepsize, constepsize;
2354 gmx_bool bDone, bAbort, do_x, do_f;
2359 int steps_accepted = 0;
2363 s_min = init_em_state();
2364 s_try = init_em_state();
2366 /* Init em and store the local state in s_try */
2367 init_em(fplog, SD, cr, inputrec,
2368 state_global, top_global, s_try, &top, &f, &f_global,
2369 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2370 nfile, fnm, &outf, &mdebin);
2372 /* Print to log file */
2373 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2375 /* Set variables for stepsize (in nm). This is the largest
2376 * step that we are going to make in any direction.
2378 ustep = inputrec->em_stepsize;
2381 /* Max number of steps */
2382 nsteps = inputrec->nsteps;
2386 /* Print to the screen */
2387 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2391 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2394 /**** HERE STARTS THE LOOP ****
2395 * count is the counter for the number of steps
2396 * bDone will be TRUE when the minimization has converged
2397 * bAbort will be TRUE when nsteps steps have been performed or when
2398 * the stepsize becomes smaller than is reasonable for machine precision
2403 while (!bDone && !bAbort)
2405 bAbort = (nsteps >= 0) && (count == nsteps);
2407 /* set new coordinates, except for first step */
2410 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2411 s_min, stepsize, s_min->f, s_try,
2412 constr, top, nrnb, wcycle, count);
2415 evaluate_energy(fplog, cr,
2416 top_global, s_try, top,
2417 inputrec, nrnb, wcycle, gstat,
2418 vsite, constr, fcd, graph, mdatoms, fr,
2419 mu_tot, enerd, vir, pres, count, count == 0);
2423 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2428 s_min->epot = s_try->epot + 1;
2431 /* Print it if necessary */
2436 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2437 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2438 (s_try->epot < s_min->epot) ? '\n' : '\r');
2441 if (s_try->epot < s_min->epot)
2443 /* Store the new (lower) energies */
2444 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2445 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2446 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2447 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2448 do_per_step(steps_accepted, inputrec->nstdisreout),
2449 do_per_step(steps_accepted, inputrec->nstorireout),
2450 fplog, count, count, eprNORMAL, TRUE,
2451 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2456 /* Now if the new energy is smaller than the previous...
2457 * or if this is the first step!
2458 * or if we did random steps!
2461 if ( (count == 0) || (s_try->epot < s_min->epot) )
2465 /* Test whether the convergence criterion is met... */
2466 bDone = (s_try->fmax < inputrec->em_tol);
2468 /* Copy the arrays for force, positions and energy */
2469 /* The 'Min' array always holds the coords and forces of the minimal
2471 swap_em_state(s_min, s_try);
2477 /* Write to trn, if necessary */
2478 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2479 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2480 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2481 top_global, inputrec, count,
2482 s_min, state_global, f_global);
2486 /* If energy is not smaller make the step smaller... */
2489 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2491 /* Reload the old state */
2492 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2493 s_min, top, mdatoms, fr, vsite, constr,
2498 /* Determine new step */
2499 stepsize = ustep/s_min->fmax;
2501 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2503 if (count == nsteps || ustep < 1e-12)
2505 if (count == nsteps || ustep < 1e-6)
2510 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2511 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2517 } /* End of the loop */
2519 /* Print some shit... */
2522 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2524 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2525 top_global, inputrec, count,
2526 s_min, state_global, f_global);
2528 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2532 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2533 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2534 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2535 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2538 finish_em(cr, outf, walltime_accounting, wcycle);
2540 /* To print the actual number of steps we needed somewhere */
2541 inputrec->nsteps = count;
2543 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2546 } /* That's all folks */
2549 double do_nm(FILE *fplog, t_commrec *cr,
2550 int nfile, const t_filenm fnm[],
2551 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2552 int gmx_unused nstglobalcomm,
2553 gmx_vsite_t *vsite, gmx_constr_t constr,
2554 int gmx_unused stepout,
2555 t_inputrec *inputrec,
2556 gmx_mtop_t *top_global, t_fcdata *fcd,
2557 t_state *state_global,
2559 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2560 gmx_edsam_t gmx_unused ed,
2562 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2563 gmx_membed_t gmx_unused membed,
2564 real gmx_unused cpt_period, real gmx_unused max_hours,
2565 const char gmx_unused *deviceOptions,
2566 unsigned long gmx_unused Flags,
2567 gmx_walltime_accounting_t walltime_accounting)
2569 const char *NM = "Normal Mode Analysis";
2571 int natoms, atom, d;
2574 gmx_localtop_t *top;
2575 gmx_enerdata_t *enerd;
2577 gmx_global_stat_t gstat;
2579 real t, t0, lambda, lam0;
2584 gmx_bool bSparse; /* use sparse matrix storage format */
2586 gmx_sparsematrix_t * sparse_matrix = NULL;
2587 real * full_matrix = NULL;
2588 em_state_t * state_work;
2590 /* added with respect to mdrun */
2591 int i, j, k, row, col;
2592 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2598 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2601 state_work = init_em_state();
2603 /* Init em and store the local state in state_minimum */
2604 init_em(fplog, NM, cr, inputrec,
2605 state_global, top_global, state_work, &top,
2607 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2608 nfile, fnm, &outf, NULL);
2610 natoms = top_global->natoms;
2618 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2619 " which MIGHT not be accurate enough for normal mode analysis.\n"
2620 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2621 " are fairly modest even if you recompile in double precision.\n\n");
2625 /* Check if we can/should use sparse storage format.
2627 * Sparse format is only useful when the Hessian itself is sparse, which it
2628 * will be when we use a cutoff.
2629 * For small systems (n<1000) it is easier to always use full matrix format, though.
2631 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2633 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2636 else if (top_global->natoms < 1000)
2638 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2643 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2649 sz = DIM*top_global->natoms;
2651 fprintf(stderr, "Allocating Hessian memory...\n\n");
2655 sparse_matrix = gmx_sparsematrix_init(sz);
2656 sparse_matrix->compressed_symmetric = TRUE;
2660 snew(full_matrix, sz*sz);
2664 /* Initial values */
2665 t0 = inputrec->init_t;
2666 lam0 = inputrec->fepvals->init_lambda;
2674 /* Write start time and temperature */
2675 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2677 /* fudge nr of steps to nr of atoms */
2678 inputrec->nsteps = natoms*2;
2682 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2683 *(top_global->name), (int)inputrec->nsteps);
2686 nnodes = cr->nnodes;
2688 /* Make evaluate_energy do a single node force calculation */
2690 evaluate_energy(fplog, cr,
2691 top_global, state_work, top,
2692 inputrec, nrnb, wcycle, gstat,
2693 vsite, constr, fcd, graph, mdatoms, fr,
2694 mu_tot, enerd, vir, pres, -1, TRUE);
2695 cr->nnodes = nnodes;
2697 /* if forces are not small, warn user */
2698 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2700 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2701 if (state_work->fmax > 1.0e-3)
2703 md_print_info(cr, fplog,
2704 "The force is probably not small enough to "
2705 "ensure that you are at a minimum.\n"
2706 "Be aware that negative eigenvalues may occur\n"
2707 "when the resulting matrix is diagonalized.\n\n");
2710 /***********************************************************
2712 * Loop over all pairs in matrix
2714 * do_force called twice. Once with positive and
2715 * once with negative displacement
2717 ************************************************************/
2719 /* Steps are divided one by one over the nodes */
2720 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2723 for (d = 0; d < DIM; d++)
2725 x_min = state_work->s.x[atom][d];
2727 state_work->s.x[atom][d] = x_min - der_range;
2729 /* Make evaluate_energy do a single node force calculation */
2731 evaluate_energy(fplog, cr,
2732 top_global, state_work, top,
2733 inputrec, nrnb, wcycle, gstat,
2734 vsite, constr, fcd, graph, mdatoms, fr,
2735 mu_tot, enerd, vir, pres, atom*2, FALSE);
2737 for (i = 0; i < natoms; i++)
2739 copy_rvec(state_work->f[i], fneg[i]);
2742 state_work->s.x[atom][d] = x_min + der_range;
2744 evaluate_energy(fplog, cr,
2745 top_global, state_work, top,
2746 inputrec, nrnb, wcycle, gstat,
2747 vsite, constr, fcd, graph, mdatoms, fr,
2748 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2749 cr->nnodes = nnodes;
2751 /* x is restored to original */
2752 state_work->s.x[atom][d] = x_min;
2754 for (j = 0; j < natoms; j++)
2756 for (k = 0; (k < DIM); k++)
2759 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2767 #define mpi_type MPI_DOUBLE
2769 #define mpi_type MPI_FLOAT
2771 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2772 cr->mpi_comm_mygroup);
2777 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2783 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2784 cr->mpi_comm_mygroup, &stat);
2789 row = (atom + node)*DIM + d;
2791 for (j = 0; j < natoms; j++)
2793 for (k = 0; k < DIM; k++)
2799 if (col >= row && dfdx[j][k] != 0.0)
2801 gmx_sparsematrix_increment_value(sparse_matrix,
2802 row, col, dfdx[j][k]);
2807 full_matrix[row*sz+col] = dfdx[j][k];
2814 if (bVerbose && fplog)
2819 /* write progress */
2820 if (MASTER(cr) && bVerbose)
2822 fprintf(stderr, "\rFinished step %d out of %d",
2823 min(atom+nnodes, natoms), natoms);
2830 fprintf(stderr, "\n\nWriting Hessian...\n");
2831 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2834 finish_em(cr, outf, walltime_accounting, wcycle);
2836 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);