1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
53 #include "gmx_fatal.h"
64 #include "md_support.h"
69 #include "mtop_util.h"
72 #include "gmx_omp_nthreads.h"
73 #include "md_logging.h"
75 #include "gromacs/fileio/confio.h"
76 #include "gromacs/fileio/trajectory_writing.h"
77 #include "gromacs/linearalgebra/mtxio.h"
78 #include "gromacs/linearalgebra/sparsematrix.h"
79 #include "gromacs/timing/wallcycle.h"
80 #include "gromacs/timing/walltime_accounting.h"
91 static em_state_t *init_em_state()
97 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
98 snew(ems->s.lambda, efptNR);
103 static void print_em_start(FILE *fplog,
105 gmx_walltime_accounting_t walltime_accounting,
106 gmx_wallcycle_t wcycle,
111 walltime_accounting_start(walltime_accounting);
113 sprintf(buf, "Started %s", name);
114 print_date_and_time(fplog, cr->nodeid, buf, NULL);
116 wallcycle_start(wcycle, ewcRUN);
118 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
119 gmx_wallcycle_t wcycle)
121 wallcycle_stop(wcycle, ewcRUN);
123 walltime_accounting_end(walltime_accounting);
126 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
129 fprintf(out, "%s:\n", minimizer);
130 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
131 fprintf(out, " Number of steps = %12d\n", nsteps);
134 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
140 "\nEnergy minimization reached the maximum number "
141 "of steps before the forces reached the requested "
142 "precision Fmax < %g.\n", ftol);
147 "\nEnergy minimization has stopped, but the forces have "
148 "not converged to the requested precision Fmax < %g (which "
149 "may not be possible for your system). It stopped "
150 "because the algorithm tried to make a new step whose size "
151 "was too small, or there was no change in the energy since "
152 "last step. Either way, we regard the minimization as "
153 "converged to within the available machine precision, "
154 "given your starting configuration and EM parameters.\n%s%s",
156 sizeof(real) < sizeof(double) ?
157 "\nDouble precision normally gives you higher accuracy, but "
158 "this is often not needed for preparing to run molecular "
162 "You might need to increase your constraint accuracy, or turn\n"
163 "off constraints altogether (set constraints = none in mdp file)\n" :
166 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
171 static void print_converged(FILE *fp, const char *alg, real ftol,
172 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
173 real epot, real fmax, int nfmax, real fnorm)
175 char buf[STEPSTRSIZE];
179 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
180 alg, ftol, gmx_step_str(count, buf));
182 else if (count < nsteps)
184 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
185 "but did not reach the requested Fmax < %g.\n",
186 alg, gmx_step_str(count, buf), ftol);
190 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
191 alg, ftol, gmx_step_str(count, buf));
195 fprintf(fp, "Potential Energy = %21.14e\n", epot);
196 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
197 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
199 fprintf(fp, "Potential Energy = %14.7e\n", epot);
200 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
201 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
205 static void get_f_norm_max(t_commrec *cr,
206 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
207 real *fnorm, real *fmax, int *a_fmax)
210 real fmax2, fmax2_0, fam;
211 int la_max, a_max, start, end, i, m, gf;
213 /* This routine finds the largest force and returns it.
214 * On parallel machines the global max is taken.
220 start = mdatoms->start;
221 end = mdatoms->homenr + start;
222 if (mdatoms->cFREEZE)
224 for (i = start; i < end; i++)
226 gf = mdatoms->cFREEZE[i];
228 for (m = 0; m < DIM; m++)
230 if (!opts->nFreeze[gf][m])
245 for (i = start; i < end; i++)
257 if (la_max >= 0 && DOMAINDECOMP(cr))
259 a_max = cr->dd->gatindex[la_max];
267 snew(sum, 2*cr->nnodes+1);
268 sum[2*cr->nodeid] = fmax2;
269 sum[2*cr->nodeid+1] = a_max;
270 sum[2*cr->nnodes] = fnorm2;
271 gmx_sumd(2*cr->nnodes+1, sum, cr);
272 fnorm2 = sum[2*cr->nnodes];
273 /* Determine the global maximum */
274 for (i = 0; i < cr->nnodes; i++)
276 if (sum[2*i] > fmax2)
279 a_max = (int)(sum[2*i+1] + 0.5);
287 *fnorm = sqrt(fnorm2);
299 static void get_state_f_norm_max(t_commrec *cr,
300 t_grpopts *opts, t_mdatoms *mdatoms,
303 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
306 void init_em(FILE *fplog, const char *title,
307 t_commrec *cr, t_inputrec *ir,
308 t_state *state_global, gmx_mtop_t *top_global,
309 em_state_t *ems, gmx_localtop_t **top,
310 rvec **f, rvec **f_global,
311 t_nrnb *nrnb, rvec mu_tot,
312 t_forcerec *fr, gmx_enerdata_t **enerd,
313 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
314 gmx_vsite_t *vsite, gmx_constr_t constr,
315 int nfile, const t_filenm fnm[],
316 gmx_mdoutf_t **outf, t_mdebin **mdebin)
318 int start, homenr, i;
323 fprintf(fplog, "Initiating %s\n", title);
326 state_global->ngtc = 0;
328 /* Initialize lambda variables */
329 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
333 if (DOMAINDECOMP(cr))
335 *top = dd_init_local_top(top_global);
337 dd_init_local_state(cr->dd, state_global, &ems->s);
341 /* Distribute the charge groups over the nodes from the master node */
342 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
343 state_global, top_global, ir,
344 &ems->s, &ems->f, mdatoms, *top,
345 fr, vsite, NULL, constr,
347 dd_store_state(cr->dd, &ems->s);
351 snew(*f_global, top_global->natoms);
361 snew(*f, top_global->natoms);
363 /* Just copy the state */
364 ems->s = *state_global;
365 snew(ems->s.x, ems->s.nalloc);
366 snew(ems->f, ems->s.nalloc);
367 for (i = 0; i < state_global->natoms; i++)
369 copy_rvec(state_global->x[i], ems->s.x[i]);
371 copy_mat(state_global->box, ems->s.box);
373 if (PAR(cr) && ir->eI != eiNM)
375 /* Initialize the particle decomposition and split the topology */
376 *top = split_system(fplog, top_global, ir, cr);
378 pd_cg_range(cr, &fr->cg0, &fr->hcg);
382 *top = gmx_mtop_generate_local_top(top_global, ir);
386 forcerec_set_excl_load(fr, *top, cr);
388 setup_bonded_threading(fr, &(*top)->idef);
390 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
392 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
401 pd_at_range(cr, &start, &homenr);
407 homenr = top_global->natoms;
409 atoms2md(top_global, ir, 0, NULL, start, homenr, mdatoms);
410 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
414 set_vsite_top(vsite, *top, mdatoms, cr);
420 if (ir->eConstrAlg == econtSHAKE &&
421 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
423 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
424 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
427 if (!DOMAINDECOMP(cr))
429 set_constraints(constr, *top, ir, mdatoms, cr);
432 if (!ir->bContinuation)
434 /* Constrain the starting coordinates */
436 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
437 ir, NULL, cr, -1, 0, mdatoms,
438 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
439 ems->s.lambda[efptFEP], &dvdl_constr,
440 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
446 *gstat = global_stat_init(ir);
449 *outf = init_mdoutf(nfile, fnm, 0, cr, ir, NULL);
452 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
457 /* Init bin for energy stuff */
458 *mdebin = init_mdebin((*outf)->fp_ene, top_global, ir, NULL);
462 calc_shifts(ems->s.box, fr->shift_vec);
465 static void finish_em(t_commrec *cr, gmx_mdoutf_t *outf,
466 gmx_walltime_accounting_t walltime_accounting,
467 gmx_wallcycle_t wcycle)
469 if (!(cr->duty & DUTY_PME))
471 /* Tell the PME only node to finish */
472 gmx_pme_send_finish(cr);
477 em_time_end(walltime_accounting, wcycle);
480 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
489 static void copy_em_coords(em_state_t *ems, t_state *state)
493 for (i = 0; (i < state->natoms); i++)
495 copy_rvec(ems->s.x[i], state->x[i]);
499 static void write_em_traj(FILE *fplog, t_commrec *cr,
501 gmx_bool bX, gmx_bool bF, const char *confout,
502 gmx_mtop_t *top_global,
503 t_inputrec *ir, gmx_int64_t step,
505 t_state *state_global, rvec *f_global)
509 if ((bX || bF || confout != NULL) && !DOMAINDECOMP(cr))
511 copy_em_coords(state, state_global);
518 mdof_flags |= MDOF_X;
522 mdof_flags |= MDOF_F;
524 write_traj(fplog, cr, outf, mdof_flags,
525 top_global, step, (double)step,
526 &state->s, state_global, state->f, f_global, NULL, NULL);
528 if (confout != NULL && MASTER(cr))
530 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
532 /* Make molecules whole only for confout writing */
533 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
537 write_sto_conf_mtop(confout,
538 *top_global->name, top_global,
539 state_global->x, NULL, ir->ePBC, state_global->box);
543 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
545 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
546 gmx_constr_t constr, gmx_localtop_t *top,
547 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
560 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
562 gmx_incons("state mismatch in do_em_step");
565 s2->flags = s1->flags;
567 if (s2->nalloc != s1->nalloc)
569 s2->nalloc = s1->nalloc;
570 srenew(s2->x, s1->nalloc);
571 srenew(ems2->f, s1->nalloc);
572 if (s2->flags & (1<<estCGP))
574 srenew(s2->cg_p, s1->nalloc);
578 s2->natoms = s1->natoms;
579 copy_mat(s1->box, s2->box);
580 /* Copy free energy state */
581 for (i = 0; i < efptNR; i++)
583 s2->lambda[i] = s1->lambda[i];
585 copy_mat(s1->box, s2->box);
588 end = md->start + md->homenr;
593 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
598 #pragma omp for schedule(static) nowait
599 for (i = start; i < end; i++)
605 for (m = 0; m < DIM; m++)
607 if (ir->opts.nFreeze[gf][m])
613 x2[i][m] = x1[i][m] + a*f[i][m];
618 if (s2->flags & (1<<estCGP))
620 /* Copy the CG p vector */
623 #pragma omp for schedule(static) nowait
624 for (i = start; i < end; i++)
626 copy_rvec(x1[i], x2[i]);
630 if (DOMAINDECOMP(cr))
632 s2->ddp_count = s1->ddp_count;
633 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
636 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
637 srenew(s2->cg_gl, s2->cg_gl_nalloc);
640 s2->ncg_gl = s1->ncg_gl;
641 #pragma omp for schedule(static) nowait
642 for (i = 0; i < s2->ncg_gl; i++)
644 s2->cg_gl[i] = s1->cg_gl[i];
646 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
652 wallcycle_start(wcycle, ewcCONSTR);
654 constrain(NULL, TRUE, TRUE, constr, &top->idef,
655 ir, NULL, cr, count, 0, md,
656 s1->x, s2->x, NULL, bMolPBC, s2->box,
657 s2->lambda[efptBONDED], &dvdl_constr,
658 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
659 wallcycle_stop(wcycle, ewcCONSTR);
663 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
664 gmx_mtop_t *top_global, t_inputrec *ir,
665 em_state_t *ems, gmx_localtop_t *top,
666 t_mdatoms *mdatoms, t_forcerec *fr,
667 gmx_vsite_t *vsite, gmx_constr_t constr,
668 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
670 /* Repartition the domain decomposition */
671 wallcycle_start(wcycle, ewcDOMDEC);
672 dd_partition_system(fplog, step, cr, FALSE, 1,
673 NULL, top_global, ir,
675 mdatoms, top, fr, vsite, NULL, constr,
676 nrnb, wcycle, FALSE);
677 dd_store_state(cr->dd, &ems->s);
678 wallcycle_stop(wcycle, ewcDOMDEC);
681 static void evaluate_energy(FILE *fplog, t_commrec *cr,
682 gmx_mtop_t *top_global,
683 em_state_t *ems, gmx_localtop_t *top,
684 t_inputrec *inputrec,
685 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
686 gmx_global_stat_t gstat,
687 gmx_vsite_t *vsite, gmx_constr_t constr,
689 t_graph *graph, t_mdatoms *mdatoms,
690 t_forcerec *fr, rvec mu_tot,
691 gmx_enerdata_t *enerd, tensor vir, tensor pres,
692 gmx_int64_t count, gmx_bool bFirst)
697 tensor force_vir, shake_vir, ekin;
698 real dvdl_constr, prescorr, enercorr, dvdlcorr;
701 /* Set the time to the initial time, the time does not change during EM */
702 t = inputrec->init_t;
705 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
707 /* This the first state or an old state used before the last ns */
713 if (inputrec->nstlist > 0)
717 else if (inputrec->nstlist == -1)
719 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
722 gmx_sumi(1, &nabnsb, cr);
730 construct_vsites(vsite, ems->s.x, 1, NULL,
731 top->idef.iparams, top->idef.il,
732 fr->ePBC, fr->bMolPBC, graph, cr, ems->s.box);
735 if (DOMAINDECOMP(cr))
739 /* Repartition the domain decomposition */
740 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
741 ems, top, mdatoms, fr, vsite, constr,
746 /* Calc force & energy on new trial position */
747 /* do_force always puts the charge groups in the box and shifts again
748 * We do not unshift, so molecules are always whole in congrad.c
750 do_force(fplog, cr, inputrec,
751 count, nrnb, wcycle, top, &top_global->groups,
752 ems->s.box, ems->s.x, &ems->s.hist,
753 ems->f, force_vir, mdatoms, enerd, fcd,
754 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
755 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
756 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
757 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
759 /* Clear the unused shake virial and pressure */
760 clear_mat(shake_vir);
763 /* Communicate stuff when parallel */
764 if (PAR(cr) && inputrec->eI != eiNM)
766 wallcycle_start(wcycle, ewcMoveE);
768 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
769 inputrec, NULL, NULL, NULL, 1, &terminate,
770 top_global, &ems->s, FALSE,
776 wallcycle_stop(wcycle, ewcMoveE);
779 /* Calculate long range corrections to pressure and energy */
780 calc_dispcorr(fplog, inputrec, fr, count, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
781 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
782 enerd->term[F_DISPCORR] = enercorr;
783 enerd->term[F_EPOT] += enercorr;
784 enerd->term[F_PRES] += prescorr;
785 enerd->term[F_DVDL] += dvdlcorr;
787 ems->epot = enerd->term[F_EPOT];
791 /* Project out the constraint components of the force */
792 wallcycle_start(wcycle, ewcCONSTR);
794 constrain(NULL, FALSE, FALSE, constr, &top->idef,
795 inputrec, NULL, cr, count, 0, mdatoms,
796 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
797 ems->s.lambda[efptBONDED], &dvdl_constr,
798 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
799 if (fr->bSepDVDL && fplog)
801 gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
803 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
804 m_add(force_vir, shake_vir, vir);
805 wallcycle_stop(wcycle, ewcCONSTR);
809 copy_mat(force_vir, vir);
813 enerd->term[F_PRES] =
814 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
816 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
818 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
820 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
824 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
826 em_state_t *s_min, em_state_t *s_b)
830 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
832 unsigned char *grpnrFREEZE;
836 fprintf(debug, "Doing reorder_partsum\n");
842 cgs_gl = dd_charge_groups_global(cr->dd);
843 index = cgs_gl->index;
845 /* Collect fm in a global vector fmg.
846 * This conflicts with the spirit of domain decomposition,
847 * but to fully optimize this a much more complicated algorithm is required.
849 snew(fmg, mtop->natoms);
851 ncg = s_min->s.ncg_gl;
852 cg_gl = s_min->s.cg_gl;
854 for (c = 0; c < ncg; c++)
859 for (a = a0; a < a1; a++)
861 copy_rvec(fm[i], fmg[a]);
865 gmx_sum(mtop->natoms*3, fmg[0], cr);
867 /* Now we will determine the part of the sum for the cgs in state s_b */
869 cg_gl = s_b->s.cg_gl;
873 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
874 for (c = 0; c < ncg; c++)
879 for (a = a0; a < a1; a++)
881 if (mdatoms->cFREEZE && grpnrFREEZE)
885 for (m = 0; m < DIM; m++)
887 if (!opts->nFreeze[gf][m])
889 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
901 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
903 em_state_t *s_min, em_state_t *s_b)
909 /* This is just the classical Polak-Ribiere calculation of beta;
910 * it looks a bit complicated since we take freeze groups into account,
911 * and might have to sum it in parallel runs.
914 if (!DOMAINDECOMP(cr) ||
915 (s_min->s.ddp_count == cr->dd->ddp_count &&
916 s_b->s.ddp_count == cr->dd->ddp_count))
922 /* This part of code can be incorrect with DD,
923 * since the atom ordering in s_b and s_min might differ.
925 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
927 if (mdatoms->cFREEZE)
929 gf = mdatoms->cFREEZE[i];
931 for (m = 0; m < DIM; m++)
933 if (!opts->nFreeze[gf][m])
935 sum += (fb[i][m] - fm[i][m])*fb[i][m];
942 /* We need to reorder cgs while summing */
943 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
947 gmx_sumd(1, &sum, cr);
950 return sum/sqr(s_min->fnorm);
953 double do_cg(FILE *fplog, t_commrec *cr,
954 int nfile, const t_filenm fnm[],
955 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
956 int gmx_unused nstglobalcomm,
957 gmx_vsite_t *vsite, gmx_constr_t constr,
958 int gmx_unused stepout,
959 t_inputrec *inputrec,
960 gmx_mtop_t *top_global, t_fcdata *fcd,
961 t_state *state_global,
963 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
964 gmx_edsam_t gmx_unused ed,
966 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
967 gmx_membed_t gmx_unused membed,
968 real gmx_unused cpt_period, real gmx_unused max_hours,
969 const char gmx_unused *deviceOptions,
970 unsigned long gmx_unused Flags,
971 gmx_walltime_accounting_t walltime_accounting)
973 const char *CG = "Polak-Ribiere Conjugate Gradients";
975 em_state_t *s_min, *s_a, *s_b, *s_c;
977 gmx_enerdata_t *enerd;
979 gmx_global_stat_t gstat;
981 rvec *f_global, *p, *sf, *sfm;
982 double gpa, gpb, gpc, tmp, sum[2], minstep;
985 real a, b, c, beta = 0.0;
989 gmx_bool converged, foundlower;
991 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
993 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
995 int i, m, gf, step, nminstep;
1000 s_min = init_em_state();
1001 s_a = init_em_state();
1002 s_b = init_em_state();
1003 s_c = init_em_state();
1005 /* Init em and store the local state in s_min */
1006 init_em(fplog, CG, cr, inputrec,
1007 state_global, top_global, s_min, &top, &f, &f_global,
1008 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1009 nfile, fnm, &outf, &mdebin);
1011 /* Print to log file */
1012 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1014 /* Max number of steps */
1015 number_steps = inputrec->nsteps;
1019 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1023 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1026 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1027 /* do_force always puts the charge groups in the box and shifts again
1028 * We do not unshift, so molecules are always whole in congrad.c
1030 evaluate_energy(fplog, cr,
1031 top_global, s_min, top,
1032 inputrec, nrnb, wcycle, gstat,
1033 vsite, constr, fcd, graph, mdatoms, fr,
1034 mu_tot, enerd, vir, pres, -1, TRUE);
1039 /* Copy stuff to the energy bin for easy printing etc. */
1040 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1041 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1042 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1044 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1045 print_ebin(outf->fp_ene, TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1046 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1050 /* Estimate/guess the initial stepsize */
1051 stepsize = inputrec->em_stepsize/s_min->fnorm;
1055 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1056 s_min->fmax, s_min->a_fmax+1);
1057 fprintf(stderr, " F-Norm = %12.5e\n",
1058 s_min->fnorm/sqrt(state_global->natoms));
1059 fprintf(stderr, "\n");
1060 /* and copy to the log file too... */
1061 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1062 s_min->fmax, s_min->a_fmax+1);
1063 fprintf(fplog, " F-Norm = %12.5e\n",
1064 s_min->fnorm/sqrt(state_global->natoms));
1065 fprintf(fplog, "\n");
1067 /* Start the loop over CG steps.
1068 * Each successful step is counted, and we continue until
1069 * we either converge or reach the max number of steps.
1072 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1075 /* start taking steps in a new direction
1076 * First time we enter the routine, beta=0, and the direction is
1077 * simply the negative gradient.
1080 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1085 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1087 if (mdatoms->cFREEZE)
1089 gf = mdatoms->cFREEZE[i];
1091 for (m = 0; m < DIM; m++)
1093 if (!inputrec->opts.nFreeze[gf][m])
1095 p[i][m] = sf[i][m] + beta*p[i][m];
1096 gpa -= p[i][m]*sf[i][m];
1097 /* f is negative gradient, thus the sign */
1106 /* Sum the gradient along the line across CPUs */
1109 gmx_sumd(1, &gpa, cr);
1112 /* Calculate the norm of the search vector */
1113 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1115 /* Just in case stepsize reaches zero due to numerical precision... */
1118 stepsize = inputrec->em_stepsize/pnorm;
1122 * Double check the value of the derivative in the search direction.
1123 * If it is positive it must be due to the old information in the
1124 * CG formula, so just remove that and start over with beta=0.
1125 * This corresponds to a steepest descent step.
1130 step--; /* Don't count this step since we are restarting */
1131 continue; /* Go back to the beginning of the big for-loop */
1134 /* Calculate minimum allowed stepsize, before the average (norm)
1135 * relative change in coordinate is smaller than precision
1138 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1140 for (m = 0; m < DIM; m++)
1142 tmp = fabs(s_min->s.x[i][m]);
1151 /* Add up from all CPUs */
1154 gmx_sumd(1, &minstep, cr);
1157 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1159 if (stepsize < minstep)
1165 /* Write coordinates if necessary */
1166 do_x = do_per_step(step, inputrec->nstxout);
1167 do_f = do_per_step(step, inputrec->nstfout);
1169 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1170 top_global, inputrec, step,
1171 s_min, state_global, f_global);
1173 /* Take a step downhill.
1174 * In theory, we should minimize the function along this direction.
1175 * That is quite possible, but it turns out to take 5-10 function evaluations
1176 * for each line. However, we dont really need to find the exact minimum -
1177 * it is much better to start a new CG step in a modified direction as soon
1178 * as we are close to it. This will save a lot of energy evaluations.
1180 * In practice, we just try to take a single step.
1181 * If it worked (i.e. lowered the energy), we increase the stepsize but
1182 * the continue straight to the next CG step without trying to find any minimum.
1183 * If it didn't work (higher energy), there must be a minimum somewhere between
1184 * the old position and the new one.
1186 * Due to the finite numerical accuracy, it turns out that it is a good idea
1187 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1188 * This leads to lower final energies in the tests I've done. / Erik
1190 s_a->epot = s_min->epot;
1192 c = a + stepsize; /* reference position along line is zero */
1194 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1196 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1197 s_min, top, mdatoms, fr, vsite, constr,
1201 /* Take a trial step (new coords in s_c) */
1202 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1203 constr, top, nrnb, wcycle, -1);
1206 /* Calculate energy for the trial step */
1207 evaluate_energy(fplog, cr,
1208 top_global, s_c, top,
1209 inputrec, nrnb, wcycle, gstat,
1210 vsite, constr, fcd, graph, mdatoms, fr,
1211 mu_tot, enerd, vir, pres, -1, FALSE);
1213 /* Calc derivative along line */
1217 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1219 for (m = 0; m < DIM; m++)
1221 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1224 /* Sum the gradient along the line across CPUs */
1227 gmx_sumd(1, &gpc, cr);
1230 /* This is the max amount of increase in energy we tolerate */
1231 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1233 /* Accept the step if the energy is lower, or if it is not significantly higher
1234 * and the line derivative is still negative.
1236 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1239 /* Great, we found a better energy. Increase step for next iteration
1240 * if we are still going down, decrease it otherwise
1244 stepsize *= 1.618034; /* The golden section */
1248 stepsize *= 0.618034; /* 1/golden section */
1253 /* New energy is the same or higher. We will have to do some work
1254 * to find a smaller value in the interval. Take smaller step next time!
1257 stepsize *= 0.618034;
1263 /* OK, if we didn't find a lower value we will have to locate one now - there must
1264 * be one in the interval [a=0,c].
1265 * The same thing is valid here, though: Don't spend dozens of iterations to find
1266 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1267 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1269 * I also have a safeguard for potentially really patological functions so we never
1270 * take more than 20 steps before we give up ...
1272 * If we already found a lower value we just skip this step and continue to the update.
1280 /* Select a new trial point.
1281 * If the derivatives at points a & c have different sign we interpolate to zero,
1282 * otherwise just do a bisection.
1284 if (gpa < 0 && gpc > 0)
1286 b = a + gpa*(a-c)/(gpc-gpa);
1293 /* safeguard if interpolation close to machine accuracy causes errors:
1294 * never go outside the interval
1296 if (b <= a || b >= c)
1301 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1303 /* Reload the old state */
1304 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1305 s_min, top, mdatoms, fr, vsite, constr,
1309 /* Take a trial step to this new point - new coords in s_b */
1310 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1311 constr, top, nrnb, wcycle, -1);
1314 /* Calculate energy for the trial step */
1315 evaluate_energy(fplog, cr,
1316 top_global, s_b, top,
1317 inputrec, nrnb, wcycle, gstat,
1318 vsite, constr, fcd, graph, mdatoms, fr,
1319 mu_tot, enerd, vir, pres, -1, FALSE);
1321 /* p does not change within a step, but since the domain decomposition
1322 * might change, we have to use cg_p of s_b here.
1327 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1329 for (m = 0; m < DIM; m++)
1331 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1334 /* Sum the gradient along the line across CPUs */
1337 gmx_sumd(1, &gpb, cr);
1342 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1343 s_a->epot, s_b->epot, s_c->epot, gpb);
1346 epot_repl = s_b->epot;
1348 /* Keep one of the intervals based on the value of the derivative at the new point */
1351 /* Replace c endpoint with b */
1352 swap_em_state(s_b, s_c);
1358 /* Replace a endpoint with b */
1359 swap_em_state(s_b, s_a);
1365 * Stop search as soon as we find a value smaller than the endpoints.
1366 * Never run more than 20 steps, no matter what.
1370 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1373 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1376 /* OK. We couldn't find a significantly lower energy.
1377 * If beta==0 this was steepest descent, and then we give up.
1378 * If not, set beta=0 and restart with steepest descent before quitting.
1388 /* Reset memory before giving up */
1394 /* Select min energy state of A & C, put the best in B.
1396 if (s_c->epot < s_a->epot)
1400 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1401 s_c->epot, s_a->epot);
1403 swap_em_state(s_b, s_c);
1411 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1412 s_a->epot, s_c->epot);
1414 swap_em_state(s_b, s_a);
1424 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1427 swap_em_state(s_b, s_c);
1432 /* new search direction */
1433 /* beta = 0 means forget all memory and restart with steepest descents. */
1434 if (nstcg && ((step % nstcg) == 0))
1440 /* s_min->fnorm cannot be zero, because then we would have converged
1444 /* Polak-Ribiere update.
1445 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1447 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1449 /* Limit beta to prevent oscillations */
1450 if (fabs(beta) > 5.0)
1456 /* update positions */
1457 swap_em_state(s_min, s_b);
1460 /* Print it if necessary */
1465 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1466 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1467 s_min->fmax, s_min->a_fmax+1);
1469 /* Store the new (lower) energies */
1470 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1471 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1472 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1474 do_log = do_per_step(step, inputrec->nstlog);
1475 do_ene = do_per_step(step, inputrec->nstenergy);
1478 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1480 print_ebin(outf->fp_ene, do_ene, FALSE, FALSE,
1481 do_log ? fplog : NULL, step, step, eprNORMAL,
1482 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1485 /* Stop when the maximum force lies below tolerance.
1486 * If we have reached machine precision, converged is already set to true.
1488 converged = converged || (s_min->fmax < inputrec->em_tol);
1490 } /* End of the loop */
1494 step--; /* we never took that last step in this case */
1497 if (s_min->fmax > inputrec->em_tol)
1501 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1502 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1509 /* If we printed energy and/or logfile last step (which was the last step)
1510 * we don't have to do it again, but otherwise print the final values.
1514 /* Write final value to log since we didn't do anything the last step */
1515 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1517 if (!do_ene || !do_log)
1519 /* Write final energy file entries */
1520 print_ebin(outf->fp_ene, !do_ene, FALSE, FALSE,
1521 !do_log ? fplog : NULL, step, step, eprNORMAL,
1522 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1526 /* Print some stuff... */
1529 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1533 * For accurate normal mode calculation it is imperative that we
1534 * store the last conformation into the full precision binary trajectory.
1536 * However, we should only do it if we did NOT already write this step
1537 * above (which we did if do_x or do_f was true).
1539 do_x = !do_per_step(step, inputrec->nstxout);
1540 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1542 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1543 top_global, inputrec, step,
1544 s_min, state_global, f_global);
1546 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1550 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1551 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1552 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1553 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1555 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1558 finish_em(cr, outf, walltime_accounting, wcycle);
1560 /* To print the actual number of steps we needed somewhere */
1561 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1564 } /* That's all folks */
1567 double do_lbfgs(FILE *fplog, t_commrec *cr,
1568 int nfile, const t_filenm fnm[],
1569 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1570 int gmx_unused nstglobalcomm,
1571 gmx_vsite_t *vsite, gmx_constr_t constr,
1572 int gmx_unused stepout,
1573 t_inputrec *inputrec,
1574 gmx_mtop_t *top_global, t_fcdata *fcd,
1577 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1578 gmx_edsam_t gmx_unused ed,
1580 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1581 gmx_membed_t gmx_unused membed,
1582 real gmx_unused cpt_period, real gmx_unused max_hours,
1583 const char gmx_unused *deviceOptions,
1584 unsigned long gmx_unused Flags,
1585 gmx_walltime_accounting_t walltime_accounting)
1587 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1589 gmx_localtop_t *top;
1590 gmx_enerdata_t *enerd;
1592 gmx_global_stat_t gstat;
1595 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1596 double stepsize, gpa, gpb, gpc, tmp, minstep;
1597 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1598 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1599 real a, b, c, maxdelta, delta;
1600 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1601 real dgdx, dgdg, sq, yr, beta;
1603 gmx_bool converged, first;
1606 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1608 int start, end, number_steps;
1610 int i, k, m, n, nfmax, gf, step;
1617 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1622 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).");
1625 n = 3*state->natoms;
1626 nmaxcorr = inputrec->nbfgscorr;
1628 /* Allocate memory */
1629 /* Use pointers to real so we dont have to loop over both atoms and
1630 * dimensions all the time...
1631 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1632 * that point to the same memory.
1645 snew(rho, nmaxcorr);
1646 snew(alpha, nmaxcorr);
1649 for (i = 0; i < nmaxcorr; i++)
1655 for (i = 0; i < nmaxcorr; i++)
1664 init_em(fplog, LBFGS, cr, inputrec,
1665 state, top_global, &ems, &top, &f, &f_global,
1666 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1667 nfile, fnm, &outf, &mdebin);
1668 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1669 * so we free some memory again.
1674 xx = (real *)state->x;
1677 start = mdatoms->start;
1678 end = mdatoms->homenr + start;
1680 /* Print to log file */
1681 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1683 do_log = do_ene = do_x = do_f = TRUE;
1685 /* Max number of steps */
1686 number_steps = inputrec->nsteps;
1688 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1690 for (i = start; i < end; i++)
1692 if (mdatoms->cFREEZE)
1694 gf = mdatoms->cFREEZE[i];
1696 for (m = 0; m < DIM; m++)
1698 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1703 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1707 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1712 construct_vsites(vsite, state->x, 1, NULL,
1713 top->idef.iparams, top->idef.il,
1714 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1717 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1718 /* do_force always puts the charge groups in the box and shifts again
1719 * We do not unshift, so molecules are always whole
1724 evaluate_energy(fplog, cr,
1725 top_global, &ems, top,
1726 inputrec, nrnb, wcycle, gstat,
1727 vsite, constr, fcd, graph, mdatoms, fr,
1728 mu_tot, enerd, vir, pres, -1, TRUE);
1733 /* Copy stuff to the energy bin for easy printing etc. */
1734 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1735 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1736 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1738 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1739 print_ebin(outf->fp_ene, TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1740 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1744 /* This is the starting energy */
1745 Epot = enerd->term[F_EPOT];
1751 /* Set the initial step.
1752 * since it will be multiplied by the non-normalized search direction
1753 * vector (force vector the first time), we scale it by the
1754 * norm of the force.
1759 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1760 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1761 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1762 fprintf(stderr, "\n");
1763 /* and copy to the log file too... */
1764 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1765 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1766 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1767 fprintf(fplog, "\n");
1771 for (i = 0; i < n; i++)
1775 dx[point][i] = ff[i]; /* Initial search direction */
1783 stepsize = 1.0/fnorm;
1786 /* Start the loop over BFGS steps.
1787 * Each successful step is counted, and we continue until
1788 * we either converge or reach the max number of steps.
1793 /* Set the gradient from the force */
1795 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1798 /* Write coordinates if necessary */
1799 do_x = do_per_step(step, inputrec->nstxout);
1800 do_f = do_per_step(step, inputrec->nstfout);
1805 mdof_flags |= MDOF_X;
1810 mdof_flags |= MDOF_F;
1813 write_traj(fplog, cr, outf, mdof_flags,
1814 top_global, step, (real)step, state, state, f, f, NULL, NULL);
1816 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1818 /* pointer to current direction - point=0 first time here */
1821 /* calculate line gradient */
1822 for (gpa = 0, i = 0; i < n; i++)
1827 /* Calculate minimum allowed stepsize, before the average (norm)
1828 * relative change in coordinate is smaller than precision
1830 for (minstep = 0, i = 0; i < n; i++)
1840 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1842 if (stepsize < minstep)
1848 /* Store old forces and coordinates */
1849 for (i = 0; i < n; i++)
1858 for (i = 0; i < n; i++)
1863 /* Take a step downhill.
1864 * In theory, we should minimize the function along this direction.
1865 * That is quite possible, but it turns out to take 5-10 function evaluations
1866 * for each line. However, we dont really need to find the exact minimum -
1867 * it is much better to start a new BFGS step in a modified direction as soon
1868 * as we are close to it. This will save a lot of energy evaluations.
1870 * In practice, we just try to take a single step.
1871 * If it worked (i.e. lowered the energy), we increase the stepsize but
1872 * the continue straight to the next BFGS step without trying to find any minimum.
1873 * If it didn't work (higher energy), there must be a minimum somewhere between
1874 * the old position and the new one.
1876 * Due to the finite numerical accuracy, it turns out that it is a good idea
1877 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1878 * This leads to lower final energies in the tests I've done. / Erik
1883 c = a + stepsize; /* reference position along line is zero */
1885 /* Check stepsize first. We do not allow displacements
1886 * larger than emstep.
1892 for (i = 0; i < n; i++)
1895 if (delta > maxdelta)
1900 if (maxdelta > inputrec->em_stepsize)
1905 while (maxdelta > inputrec->em_stepsize);
1907 /* Take a trial step */
1908 for (i = 0; i < n; i++)
1910 xc[i] = lastx[i] + c*s[i];
1914 /* Calculate energy for the trial step */
1915 ems.s.x = (rvec *)xc;
1917 evaluate_energy(fplog, cr,
1918 top_global, &ems, top,
1919 inputrec, nrnb, wcycle, gstat,
1920 vsite, constr, fcd, graph, mdatoms, fr,
1921 mu_tot, enerd, vir, pres, step, FALSE);
1924 /* Calc derivative along line */
1925 for (gpc = 0, i = 0; i < n; i++)
1927 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1929 /* Sum the gradient along the line across CPUs */
1932 gmx_sumd(1, &gpc, cr);
1935 /* This is the max amount of increase in energy we tolerate */
1936 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1938 /* Accept the step if the energy is lower, or if it is not significantly higher
1939 * and the line derivative is still negative.
1941 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1944 /* Great, we found a better energy. Increase step for next iteration
1945 * if we are still going down, decrease it otherwise
1949 stepsize *= 1.618034; /* The golden section */
1953 stepsize *= 0.618034; /* 1/golden section */
1958 /* New energy is the same or higher. We will have to do some work
1959 * to find a smaller value in the interval. Take smaller step next time!
1962 stepsize *= 0.618034;
1965 /* OK, if we didn't find a lower value we will have to locate one now - there must
1966 * be one in the interval [a=0,c].
1967 * The same thing is valid here, though: Don't spend dozens of iterations to find
1968 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1969 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1971 * I also have a safeguard for potentially really patological functions so we never
1972 * take more than 20 steps before we give up ...
1974 * If we already found a lower value we just skip this step and continue to the update.
1983 /* Select a new trial point.
1984 * If the derivatives at points a & c have different sign we interpolate to zero,
1985 * otherwise just do a bisection.
1988 if (gpa < 0 && gpc > 0)
1990 b = a + gpa*(a-c)/(gpc-gpa);
1997 /* safeguard if interpolation close to machine accuracy causes errors:
1998 * never go outside the interval
2000 if (b <= a || b >= c)
2005 /* Take a trial step */
2006 for (i = 0; i < n; i++)
2008 xb[i] = lastx[i] + b*s[i];
2012 /* Calculate energy for the trial step */
2013 ems.s.x = (rvec *)xb;
2015 evaluate_energy(fplog, cr,
2016 top_global, &ems, top,
2017 inputrec, nrnb, wcycle, gstat,
2018 vsite, constr, fcd, graph, mdatoms, fr,
2019 mu_tot, enerd, vir, pres, step, FALSE);
2024 for (gpb = 0, i = 0; i < n; i++)
2026 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2029 /* Sum the gradient along the line across CPUs */
2032 gmx_sumd(1, &gpb, cr);
2035 /* Keep one of the intervals based on the value of the derivative at the new point */
2038 /* Replace c endpoint with b */
2042 /* swap coord pointers b/c */
2052 /* Replace a endpoint with b */
2056 /* swap coord pointers a/b */
2066 * Stop search as soon as we find a value smaller than the endpoints,
2067 * or if the tolerance is below machine precision.
2068 * Never run more than 20 steps, no matter what.
2072 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2074 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2076 /* OK. We couldn't find a significantly lower energy.
2077 * If ncorr==0 this was steepest descent, and then we give up.
2078 * If not, reset memory to restart as steepest descent before quitting.
2090 /* Search in gradient direction */
2091 for (i = 0; i < n; i++)
2093 dx[point][i] = ff[i];
2095 /* Reset stepsize */
2096 stepsize = 1.0/fnorm;
2101 /* Select min energy state of A & C, put the best in xx/ff/Epot
2107 for (i = 0; i < n; i++)
2118 for (i = 0; i < n; i++)
2132 for (i = 0; i < n; i++)
2140 /* Update the memory information, and calculate a new
2141 * approximation of the inverse hessian
2144 /* Have new data in Epot, xx, ff */
2145 if (ncorr < nmaxcorr)
2150 for (i = 0; i < n; i++)
2152 dg[point][i] = lastf[i]-ff[i];
2153 dx[point][i] *= stepsize;
2158 for (i = 0; i < n; i++)
2160 dgdg += dg[point][i]*dg[point][i];
2161 dgdx += dg[point][i]*dx[point][i];
2166 rho[point] = 1.0/dgdx;
2169 if (point >= nmaxcorr)
2175 for (i = 0; i < n; i++)
2182 /* Recursive update. First go back over the memory points */
2183 for (k = 0; k < ncorr; k++)
2192 for (i = 0; i < n; i++)
2194 sq += dx[cp][i]*p[i];
2197 alpha[cp] = rho[cp]*sq;
2199 for (i = 0; i < n; i++)
2201 p[i] -= alpha[cp]*dg[cp][i];
2205 for (i = 0; i < n; i++)
2210 /* And then go forward again */
2211 for (k = 0; k < ncorr; k++)
2214 for (i = 0; i < n; i++)
2216 yr += p[i]*dg[cp][i];
2220 beta = alpha[cp]-beta;
2222 for (i = 0; i < n; i++)
2224 p[i] += beta*dx[cp][i];
2234 for (i = 0; i < n; i++)
2238 dx[point][i] = p[i];
2248 /* Test whether the convergence criterion is met */
2249 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2251 /* Print it if necessary */
2256 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2257 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2259 /* Store the new (lower) energies */
2260 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2261 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2262 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2263 do_log = do_per_step(step, inputrec->nstlog);
2264 do_ene = do_per_step(step, inputrec->nstenergy);
2267 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2269 print_ebin(outf->fp_ene, do_ene, FALSE, FALSE,
2270 do_log ? fplog : NULL, step, step, eprNORMAL,
2271 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2274 /* Stop when the maximum force lies below tolerance.
2275 * If we have reached machine precision, converged is already set to true.
2278 converged = converged || (fmax < inputrec->em_tol);
2280 } /* End of the loop */
2284 step--; /* we never took that last step in this case */
2287 if (fmax > inputrec->em_tol)
2291 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2292 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2297 /* If we printed energy and/or logfile last step (which was the last step)
2298 * we don't have to do it again, but otherwise print the final values.
2300 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2302 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2304 if (!do_ene || !do_log) /* Write final energy file entries */
2306 print_ebin(outf->fp_ene, !do_ene, FALSE, FALSE,
2307 !do_log ? fplog : NULL, step, step, eprNORMAL,
2308 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2311 /* Print some stuff... */
2314 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2318 * For accurate normal mode calculation it is imperative that we
2319 * store the last conformation into the full precision binary trajectory.
2321 * However, we should only do it if we did NOT already write this step
2322 * above (which we did if do_x or do_f was true).
2324 do_x = !do_per_step(step, inputrec->nstxout);
2325 do_f = !do_per_step(step, inputrec->nstfout);
2326 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2327 top_global, inputrec, step,
2332 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2333 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2334 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2335 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2337 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2340 finish_em(cr, outf, walltime_accounting, wcycle);
2342 /* To print the actual number of steps we needed somewhere */
2343 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2346 } /* That's all folks */
2349 double do_steep(FILE *fplog, t_commrec *cr,
2350 int nfile, const t_filenm fnm[],
2351 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2352 int gmx_unused nstglobalcomm,
2353 gmx_vsite_t *vsite, gmx_constr_t constr,
2354 int gmx_unused stepout,
2355 t_inputrec *inputrec,
2356 gmx_mtop_t *top_global, t_fcdata *fcd,
2357 t_state *state_global,
2359 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2360 gmx_edsam_t gmx_unused ed,
2362 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2363 gmx_membed_t gmx_unused membed,
2364 real gmx_unused cpt_period, real gmx_unused max_hours,
2365 const char gmx_unused *deviceOptions,
2366 unsigned long gmx_unused Flags,
2367 gmx_walltime_accounting_t walltime_accounting)
2369 const char *SD = "Steepest Descents";
2370 em_state_t *s_min, *s_try;
2372 gmx_localtop_t *top;
2373 gmx_enerdata_t *enerd;
2375 gmx_global_stat_t gstat;
2377 real stepsize, constepsize;
2381 gmx_bool bDone, bAbort, do_x, do_f;
2386 int steps_accepted = 0;
2390 s_min = init_em_state();
2391 s_try = init_em_state();
2393 /* Init em and store the local state in s_try */
2394 init_em(fplog, SD, cr, inputrec,
2395 state_global, top_global, s_try, &top, &f, &f_global,
2396 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2397 nfile, fnm, &outf, &mdebin);
2399 /* Print to log file */
2400 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2402 /* Set variables for stepsize (in nm). This is the largest
2403 * step that we are going to make in any direction.
2405 ustep = inputrec->em_stepsize;
2408 /* Max number of steps */
2409 nsteps = inputrec->nsteps;
2413 /* Print to the screen */
2414 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2418 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2421 /**** HERE STARTS THE LOOP ****
2422 * count is the counter for the number of steps
2423 * bDone will be TRUE when the minimization has converged
2424 * bAbort will be TRUE when nsteps steps have been performed or when
2425 * the stepsize becomes smaller than is reasonable for machine precision
2430 while (!bDone && !bAbort)
2432 bAbort = (nsteps >= 0) && (count == nsteps);
2434 /* set new coordinates, except for first step */
2437 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2438 s_min, stepsize, s_min->f, s_try,
2439 constr, top, nrnb, wcycle, count);
2442 evaluate_energy(fplog, cr,
2443 top_global, s_try, top,
2444 inputrec, nrnb, wcycle, gstat,
2445 vsite, constr, fcd, graph, mdatoms, fr,
2446 mu_tot, enerd, vir, pres, count, count == 0);
2450 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2455 s_min->epot = s_try->epot + 1;
2458 /* Print it if necessary */
2463 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2464 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2465 (s_try->epot < s_min->epot) ? '\n' : '\r');
2468 if (s_try->epot < s_min->epot)
2470 /* Store the new (lower) energies */
2471 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2472 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2473 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2474 print_ebin(outf->fp_ene, TRUE,
2475 do_per_step(steps_accepted, inputrec->nstdisreout),
2476 do_per_step(steps_accepted, inputrec->nstorireout),
2477 fplog, count, count, eprNORMAL, TRUE,
2478 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2483 /* Now if the new energy is smaller than the previous...
2484 * or if this is the first step!
2485 * or if we did random steps!
2488 if ( (count == 0) || (s_try->epot < s_min->epot) )
2492 /* Test whether the convergence criterion is met... */
2493 bDone = (s_try->fmax < inputrec->em_tol);
2495 /* Copy the arrays for force, positions and energy */
2496 /* The 'Min' array always holds the coords and forces of the minimal
2498 swap_em_state(s_min, s_try);
2504 /* Write to trn, if necessary */
2505 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2506 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2507 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2508 top_global, inputrec, count,
2509 s_min, state_global, f_global);
2513 /* If energy is not smaller make the step smaller... */
2516 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2518 /* Reload the old state */
2519 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2520 s_min, top, mdatoms, fr, vsite, constr,
2525 /* Determine new step */
2526 stepsize = ustep/s_min->fmax;
2528 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2530 if (count == nsteps || ustep < 1e-12)
2532 if (count == nsteps || ustep < 1e-6)
2537 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2538 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2544 } /* End of the loop */
2546 /* Print some shit... */
2549 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2551 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2552 top_global, inputrec, count,
2553 s_min, state_global, f_global);
2555 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2559 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2560 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2561 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2562 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2565 finish_em(cr, outf, walltime_accounting, wcycle);
2567 /* To print the actual number of steps we needed somewhere */
2568 inputrec->nsteps = count;
2570 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2573 } /* That's all folks */
2576 double do_nm(FILE *fplog, t_commrec *cr,
2577 int nfile, const t_filenm fnm[],
2578 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2579 int gmx_unused nstglobalcomm,
2580 gmx_vsite_t *vsite, gmx_constr_t constr,
2581 int gmx_unused stepout,
2582 t_inputrec *inputrec,
2583 gmx_mtop_t *top_global, t_fcdata *fcd,
2584 t_state *state_global,
2586 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2587 gmx_edsam_t gmx_unused ed,
2589 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2590 gmx_membed_t gmx_unused membed,
2591 real gmx_unused cpt_period, real gmx_unused max_hours,
2592 const char gmx_unused *deviceOptions,
2593 unsigned long gmx_unused Flags,
2594 gmx_walltime_accounting_t walltime_accounting)
2596 const char *NM = "Normal Mode Analysis";
2598 int natoms, atom, d;
2601 gmx_localtop_t *top;
2602 gmx_enerdata_t *enerd;
2604 gmx_global_stat_t gstat;
2606 real t, t0, lambda, lam0;
2611 gmx_bool bSparse; /* use sparse matrix storage format */
2613 gmx_sparsematrix_t * sparse_matrix = NULL;
2614 real * full_matrix = NULL;
2615 em_state_t * state_work;
2617 /* added with respect to mdrun */
2618 int i, j, k, row, col;
2619 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2625 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2628 state_work = init_em_state();
2630 /* Init em and store the local state in state_minimum */
2631 init_em(fplog, NM, cr, inputrec,
2632 state_global, top_global, state_work, &top,
2634 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2635 nfile, fnm, &outf, NULL);
2637 natoms = top_global->natoms;
2645 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2646 " which MIGHT not be accurate enough for normal mode analysis.\n"
2647 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2648 " are fairly modest even if you recompile in double precision.\n\n");
2652 /* Check if we can/should use sparse storage format.
2654 * Sparse format is only useful when the Hessian itself is sparse, which it
2655 * will be when we use a cutoff.
2656 * For small systems (n<1000) it is easier to always use full matrix format, though.
2658 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2660 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2663 else if (top_global->natoms < 1000)
2665 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2670 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2676 sz = DIM*top_global->natoms;
2678 fprintf(stderr, "Allocating Hessian memory...\n\n");
2682 sparse_matrix = gmx_sparsematrix_init(sz);
2683 sparse_matrix->compressed_symmetric = TRUE;
2687 snew(full_matrix, sz*sz);
2691 /* Initial values */
2692 t0 = inputrec->init_t;
2693 lam0 = inputrec->fepvals->init_lambda;
2701 /* Write start time and temperature */
2702 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2704 /* fudge nr of steps to nr of atoms */
2705 inputrec->nsteps = natoms*2;
2709 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2710 *(top_global->name), (int)inputrec->nsteps);
2713 nnodes = cr->nnodes;
2715 /* Make evaluate_energy do a single node force calculation */
2717 evaluate_energy(fplog, cr,
2718 top_global, state_work, top,
2719 inputrec, nrnb, wcycle, gstat,
2720 vsite, constr, fcd, graph, mdatoms, fr,
2721 mu_tot, enerd, vir, pres, -1, TRUE);
2722 cr->nnodes = nnodes;
2724 /* if forces are not small, warn user */
2725 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2727 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2728 if (state_work->fmax > 1.0e-3)
2730 md_print_info(cr, fplog,
2731 "The force is probably not small enough to "
2732 "ensure that you are at a minimum.\n"
2733 "Be aware that negative eigenvalues may occur\n"
2734 "when the resulting matrix is diagonalized.\n\n");
2737 /***********************************************************
2739 * Loop over all pairs in matrix
2741 * do_force called twice. Once with positive and
2742 * once with negative displacement
2744 ************************************************************/
2746 /* Steps are divided one by one over the nodes */
2747 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2750 for (d = 0; d < DIM; d++)
2752 x_min = state_work->s.x[atom][d];
2754 state_work->s.x[atom][d] = x_min - der_range;
2756 /* Make evaluate_energy do a single node force calculation */
2758 evaluate_energy(fplog, cr,
2759 top_global, state_work, top,
2760 inputrec, nrnb, wcycle, gstat,
2761 vsite, constr, fcd, graph, mdatoms, fr,
2762 mu_tot, enerd, vir, pres, atom*2, FALSE);
2764 for (i = 0; i < natoms; i++)
2766 copy_rvec(state_work->f[i], fneg[i]);
2769 state_work->s.x[atom][d] = x_min + der_range;
2771 evaluate_energy(fplog, cr,
2772 top_global, state_work, top,
2773 inputrec, nrnb, wcycle, gstat,
2774 vsite, constr, fcd, graph, mdatoms, fr,
2775 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2776 cr->nnodes = nnodes;
2778 /* x is restored to original */
2779 state_work->s.x[atom][d] = x_min;
2781 for (j = 0; j < natoms; j++)
2783 for (k = 0; (k < DIM); k++)
2786 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2794 #define mpi_type MPI_DOUBLE
2796 #define mpi_type MPI_FLOAT
2798 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2799 cr->mpi_comm_mygroup);
2804 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2810 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2811 cr->mpi_comm_mygroup, &stat);
2816 row = (atom + node)*DIM + d;
2818 for (j = 0; j < natoms; j++)
2820 for (k = 0; k < DIM; k++)
2826 if (col >= row && dfdx[j][k] != 0.0)
2828 gmx_sparsematrix_increment_value(sparse_matrix,
2829 row, col, dfdx[j][k]);
2834 full_matrix[row*sz+col] = dfdx[j][k];
2841 if (bVerbose && fplog)
2846 /* write progress */
2847 if (MASTER(cr) && bVerbose)
2849 fprintf(stderr, "\rFinished step %d out of %d",
2850 min(atom+nnodes, natoms), natoms);
2857 fprintf(stderr, "\n\nWriting Hessian...\n");
2858 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2861 finish_em(cr, outf, walltime_accounting, wcycle);
2863 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);