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.
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,
109 walltime_accounting_start(walltime_accounting);
110 wallcycle_start(wcycle, ewcRUN);
111 print_start(fplog, cr, walltime_accounting, name);
113 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
114 gmx_wallcycle_t wcycle)
116 wallcycle_stop(wcycle, ewcRUN);
118 walltime_accounting_end(walltime_accounting);
121 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
124 fprintf(out, "%s:\n", minimizer);
125 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
126 fprintf(out, " Number of steps = %12d\n", nsteps);
129 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
135 "\nEnergy minimization reached the maximum number "
136 "of steps before the forces reached the requested "
137 "precision Fmax < %g.\n", ftol);
142 "\nEnergy minimization has stopped, but the forces have "
143 "not converged to the requested precision Fmax < %g (which "
144 "may not be possible for your system). It stopped "
145 "because the algorithm tried to make a new step whose size "
146 "was too small, or there was no change in the energy since "
147 "last step. Either way, we regard the minimization as "
148 "converged to within the available machine precision, "
149 "given your starting configuration and EM parameters.\n%s%s",
151 sizeof(real) < sizeof(double) ?
152 "\nDouble precision normally gives you higher accuracy, but "
153 "this is often not needed for preparing to run molecular "
157 "You might need to increase your constraint accuracy, or turn\n"
158 "off constraints altogether (set constraints = none in mdp file)\n" :
161 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
166 static void print_converged(FILE *fp, const char *alg, real ftol,
167 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
168 real epot, real fmax, int nfmax, real fnorm)
170 char buf[STEPSTRSIZE];
174 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
175 alg, ftol, gmx_step_str(count, buf));
177 else if (count < nsteps)
179 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
180 "but did not reach the requested Fmax < %g.\n",
181 alg, gmx_step_str(count, buf), ftol);
185 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
186 alg, ftol, gmx_step_str(count, buf));
190 fprintf(fp, "Potential Energy = %21.14e\n", epot);
191 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
192 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
194 fprintf(fp, "Potential Energy = %14.7e\n", epot);
195 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
196 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
200 static void get_f_norm_max(t_commrec *cr,
201 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
202 real *fnorm, real *fmax, int *a_fmax)
205 real fmax2, fmax2_0, fam;
206 int la_max, a_max, start, end, i, m, gf;
208 /* This routine finds the largest force and returns it.
209 * On parallel machines the global max is taken.
216 end = mdatoms->homenr;
217 if (mdatoms->cFREEZE)
219 for (i = start; i < end; i++)
221 gf = mdatoms->cFREEZE[i];
223 for (m = 0; m < DIM; m++)
225 if (!opts->nFreeze[gf][m])
240 for (i = start; i < end; i++)
252 if (la_max >= 0 && DOMAINDECOMP(cr))
254 a_max = cr->dd->gatindex[la_max];
262 snew(sum, 2*cr->nnodes+1);
263 sum[2*cr->nodeid] = fmax2;
264 sum[2*cr->nodeid+1] = a_max;
265 sum[2*cr->nnodes] = fnorm2;
266 gmx_sumd(2*cr->nnodes+1, sum, cr);
267 fnorm2 = sum[2*cr->nnodes];
268 /* Determine the global maximum */
269 for (i = 0; i < cr->nnodes; i++)
271 if (sum[2*i] > fmax2)
274 a_max = (int)(sum[2*i+1] + 0.5);
282 *fnorm = sqrt(fnorm2);
294 static void get_state_f_norm_max(t_commrec *cr,
295 t_grpopts *opts, t_mdatoms *mdatoms,
298 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
301 void init_em(FILE *fplog, const char *title,
302 t_commrec *cr, t_inputrec *ir,
303 t_state *state_global, gmx_mtop_t *top_global,
304 em_state_t *ems, gmx_localtop_t **top,
305 rvec **f, rvec **f_global,
306 t_nrnb *nrnb, rvec mu_tot,
307 t_forcerec *fr, gmx_enerdata_t **enerd,
308 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
309 gmx_vsite_t *vsite, gmx_constr_t constr,
310 int nfile, const t_filenm fnm[],
311 gmx_mdoutf_t *outf, t_mdebin **mdebin)
318 fprintf(fplog, "Initiating %s\n", title);
321 state_global->ngtc = 0;
323 /* Initialize lambda variables */
324 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
328 if (DOMAINDECOMP(cr))
330 *top = dd_init_local_top(top_global);
332 dd_init_local_state(cr->dd, state_global, &ems->s);
336 /* Distribute the charge groups over the nodes from the master node */
337 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
338 state_global, top_global, ir,
339 &ems->s, &ems->f, mdatoms, *top,
340 fr, vsite, NULL, constr,
342 dd_store_state(cr->dd, &ems->s);
346 snew(*f_global, top_global->natoms);
356 snew(*f, top_global->natoms);
358 /* Just copy the state */
359 ems->s = *state_global;
360 snew(ems->s.x, ems->s.nalloc);
361 snew(ems->f, ems->s.nalloc);
362 for (i = 0; i < state_global->natoms; i++)
364 copy_rvec(state_global->x[i], ems->s.x[i]);
366 copy_mat(state_global->box, ems->s.box);
368 *top = gmx_mtop_generate_local_top(top_global, ir);
371 forcerec_set_excl_load(fr, *top);
373 setup_bonded_threading(fr, &(*top)->idef);
375 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
377 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
384 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
385 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
389 set_vsite_top(vsite, *top, mdatoms, cr);
395 if (ir->eConstrAlg == econtSHAKE &&
396 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
398 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
399 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
402 if (!DOMAINDECOMP(cr))
404 set_constraints(constr, *top, ir, mdatoms, cr);
407 if (!ir->bContinuation)
409 /* Constrain the starting coordinates */
411 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
412 ir, NULL, cr, -1, 0, mdatoms,
413 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
414 ems->s.lambda[efptFEP], &dvdl_constr,
415 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
421 *gstat = global_stat_init(ir);
424 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, NULL);
427 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
432 /* Init bin for energy stuff */
433 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
437 calc_shifts(ems->s.box, fr->shift_vec);
440 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
441 gmx_walltime_accounting_t walltime_accounting,
442 gmx_wallcycle_t wcycle)
444 if (!(cr->duty & DUTY_PME))
446 /* Tell the PME only node to finish */
447 gmx_pme_send_finish(cr);
452 em_time_end(walltime_accounting, wcycle);
455 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
464 static void copy_em_coords(em_state_t *ems, t_state *state)
468 for (i = 0; (i < state->natoms); i++)
470 copy_rvec(ems->s.x[i], state->x[i]);
474 static void write_em_traj(FILE *fplog, t_commrec *cr,
476 gmx_bool bX, gmx_bool bF, const char *confout,
477 gmx_mtop_t *top_global,
478 t_inputrec *ir, gmx_int64_t step,
480 t_state *state_global, rvec *f_global)
484 if ((bX || bF || confout != NULL) && !DOMAINDECOMP(cr))
486 copy_em_coords(state, state_global);
493 mdof_flags |= MDOF_X;
497 mdof_flags |= MDOF_F;
499 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
500 top_global, step, (double)step,
501 &state->s, state_global, state->f, f_global);
503 if (confout != NULL && MASTER(cr))
505 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
507 /* Make molecules whole only for confout writing */
508 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
512 write_sto_conf_mtop(confout,
513 *top_global->name, top_global,
514 state_global->x, NULL, ir->ePBC, state_global->box);
518 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
520 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
521 gmx_constr_t constr, gmx_localtop_t *top,
522 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
535 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
537 gmx_incons("state mismatch in do_em_step");
540 s2->flags = s1->flags;
542 if (s2->nalloc != s1->nalloc)
544 s2->nalloc = s1->nalloc;
545 srenew(s2->x, s1->nalloc);
546 srenew(ems2->f, s1->nalloc);
547 if (s2->flags & (1<<estCGP))
549 srenew(s2->cg_p, s1->nalloc);
553 s2->natoms = s1->natoms;
554 copy_mat(s1->box, s2->box);
555 /* Copy free energy state */
556 for (i = 0; i < efptNR; i++)
558 s2->lambda[i] = s1->lambda[i];
560 copy_mat(s1->box, s2->box);
568 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
573 #pragma omp for schedule(static) nowait
574 for (i = start; i < end; i++)
580 for (m = 0; m < DIM; m++)
582 if (ir->opts.nFreeze[gf][m])
588 x2[i][m] = x1[i][m] + a*f[i][m];
593 if (s2->flags & (1<<estCGP))
595 /* Copy the CG p vector */
598 #pragma omp for schedule(static) nowait
599 for (i = start; i < end; i++)
601 copy_rvec(x1[i], x2[i]);
605 if (DOMAINDECOMP(cr))
607 s2->ddp_count = s1->ddp_count;
608 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
611 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
612 srenew(s2->cg_gl, s2->cg_gl_nalloc);
615 s2->ncg_gl = s1->ncg_gl;
616 #pragma omp for schedule(static) nowait
617 for (i = 0; i < s2->ncg_gl; i++)
619 s2->cg_gl[i] = s1->cg_gl[i];
621 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
627 wallcycle_start(wcycle, ewcCONSTR);
629 constrain(NULL, TRUE, TRUE, constr, &top->idef,
630 ir, NULL, cr, count, 0, md,
631 s1->x, s2->x, NULL, bMolPBC, s2->box,
632 s2->lambda[efptBONDED], &dvdl_constr,
633 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
634 wallcycle_stop(wcycle, ewcCONSTR);
638 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
639 gmx_mtop_t *top_global, t_inputrec *ir,
640 em_state_t *ems, gmx_localtop_t *top,
641 t_mdatoms *mdatoms, t_forcerec *fr,
642 gmx_vsite_t *vsite, gmx_constr_t constr,
643 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
645 /* Repartition the domain decomposition */
646 wallcycle_start(wcycle, ewcDOMDEC);
647 dd_partition_system(fplog, step, cr, FALSE, 1,
648 NULL, top_global, ir,
650 mdatoms, top, fr, vsite, NULL, constr,
651 nrnb, wcycle, FALSE);
652 dd_store_state(cr->dd, &ems->s);
653 wallcycle_stop(wcycle, ewcDOMDEC);
656 static void evaluate_energy(FILE *fplog, t_commrec *cr,
657 gmx_mtop_t *top_global,
658 em_state_t *ems, gmx_localtop_t *top,
659 t_inputrec *inputrec,
660 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
661 gmx_global_stat_t gstat,
662 gmx_vsite_t *vsite, gmx_constr_t constr,
664 t_graph *graph, t_mdatoms *mdatoms,
665 t_forcerec *fr, rvec mu_tot,
666 gmx_enerdata_t *enerd, tensor vir, tensor pres,
667 gmx_int64_t count, gmx_bool bFirst)
672 tensor force_vir, shake_vir, ekin;
673 real dvdl_constr, prescorr, enercorr, dvdlcorr;
676 /* Set the time to the initial time, the time does not change during EM */
677 t = inputrec->init_t;
680 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
682 /* This is the first state or an old state used before the last ns */
688 if (inputrec->nstlist > 0)
692 else if (inputrec->nstlist == -1)
694 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
697 gmx_sumi(1, &nabnsb, cr);
705 construct_vsites(vsite, ems->s.x, 1, NULL,
706 top->idef.iparams, top->idef.il,
707 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
710 if (DOMAINDECOMP(cr) && bNS)
712 /* Repartition the domain decomposition */
713 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
714 ems, top, mdatoms, fr, vsite, constr,
718 /* Calc force & energy on new trial position */
719 /* do_force always puts the charge groups in the box and shifts again
720 * We do not unshift, so molecules are always whole in congrad.c
722 do_force(fplog, cr, inputrec,
723 count, nrnb, wcycle, top, &top_global->groups,
724 ems->s.box, ems->s.x, &ems->s.hist,
725 ems->f, force_vir, mdatoms, enerd, fcd,
726 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
727 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
728 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
729 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
731 /* Clear the unused shake virial and pressure */
732 clear_mat(shake_vir);
735 /* Communicate stuff when parallel */
736 if (PAR(cr) && inputrec->eI != eiNM)
738 wallcycle_start(wcycle, ewcMoveE);
740 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
741 inputrec, NULL, NULL, NULL, 1, &terminate,
742 top_global, &ems->s, FALSE,
748 wallcycle_stop(wcycle, ewcMoveE);
751 /* Calculate long range corrections to pressure and energy */
752 calc_dispcorr(fplog, inputrec, fr, count, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
753 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
754 enerd->term[F_DISPCORR] = enercorr;
755 enerd->term[F_EPOT] += enercorr;
756 enerd->term[F_PRES] += prescorr;
757 enerd->term[F_DVDL] += dvdlcorr;
759 ems->epot = enerd->term[F_EPOT];
763 /* Project out the constraint components of the force */
764 wallcycle_start(wcycle, ewcCONSTR);
766 constrain(NULL, FALSE, FALSE, constr, &top->idef,
767 inputrec, NULL, cr, count, 0, mdatoms,
768 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
769 ems->s.lambda[efptBONDED], &dvdl_constr,
770 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
771 if (fr->bSepDVDL && fplog)
773 gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
775 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
776 m_add(force_vir, shake_vir, vir);
777 wallcycle_stop(wcycle, ewcCONSTR);
781 copy_mat(force_vir, vir);
785 enerd->term[F_PRES] =
786 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
788 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
790 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
792 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
796 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
798 em_state_t *s_min, em_state_t *s_b)
802 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
804 unsigned char *grpnrFREEZE;
808 fprintf(debug, "Doing reorder_partsum\n");
814 cgs_gl = dd_charge_groups_global(cr->dd);
815 index = cgs_gl->index;
817 /* Collect fm in a global vector fmg.
818 * This conflicts with the spirit of domain decomposition,
819 * but to fully optimize this a much more complicated algorithm is required.
821 snew(fmg, mtop->natoms);
823 ncg = s_min->s.ncg_gl;
824 cg_gl = s_min->s.cg_gl;
826 for (c = 0; c < ncg; c++)
831 for (a = a0; a < a1; a++)
833 copy_rvec(fm[i], fmg[a]);
837 gmx_sum(mtop->natoms*3, fmg[0], cr);
839 /* Now we will determine the part of the sum for the cgs in state s_b */
841 cg_gl = s_b->s.cg_gl;
845 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
846 for (c = 0; c < ncg; c++)
851 for (a = a0; a < a1; a++)
853 if (mdatoms->cFREEZE && grpnrFREEZE)
857 for (m = 0; m < DIM; m++)
859 if (!opts->nFreeze[gf][m])
861 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
873 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
875 em_state_t *s_min, em_state_t *s_b)
881 /* This is just the classical Polak-Ribiere calculation of beta;
882 * it looks a bit complicated since we take freeze groups into account,
883 * and might have to sum it in parallel runs.
886 if (!DOMAINDECOMP(cr) ||
887 (s_min->s.ddp_count == cr->dd->ddp_count &&
888 s_b->s.ddp_count == cr->dd->ddp_count))
894 /* This part of code can be incorrect with DD,
895 * since the atom ordering in s_b and s_min might differ.
897 for (i = 0; i < mdatoms->homenr; i++)
899 if (mdatoms->cFREEZE)
901 gf = mdatoms->cFREEZE[i];
903 for (m = 0; m < DIM; m++)
905 if (!opts->nFreeze[gf][m])
907 sum += (fb[i][m] - fm[i][m])*fb[i][m];
914 /* We need to reorder cgs while summing */
915 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
919 gmx_sumd(1, &sum, cr);
922 return sum/sqr(s_min->fnorm);
925 double do_cg(FILE *fplog, t_commrec *cr,
926 int nfile, const t_filenm fnm[],
927 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
928 int gmx_unused nstglobalcomm,
929 gmx_vsite_t *vsite, gmx_constr_t constr,
930 int gmx_unused stepout,
931 t_inputrec *inputrec,
932 gmx_mtop_t *top_global, t_fcdata *fcd,
933 t_state *state_global,
935 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
936 gmx_edsam_t gmx_unused ed,
938 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
939 gmx_membed_t gmx_unused membed,
940 real gmx_unused cpt_period, real gmx_unused max_hours,
941 const char gmx_unused *deviceOptions,
942 unsigned long gmx_unused Flags,
943 gmx_walltime_accounting_t walltime_accounting)
945 const char *CG = "Polak-Ribiere Conjugate Gradients";
947 em_state_t *s_min, *s_a, *s_b, *s_c;
949 gmx_enerdata_t *enerd;
951 gmx_global_stat_t gstat;
953 rvec *f_global, *p, *sf, *sfm;
954 double gpa, gpb, gpc, tmp, sum[2], minstep;
957 real a, b, c, beta = 0.0;
961 gmx_bool converged, foundlower;
963 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
965 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
967 int i, m, gf, step, nminstep;
972 s_min = init_em_state();
973 s_a = init_em_state();
974 s_b = init_em_state();
975 s_c = init_em_state();
977 /* Init em and store the local state in s_min */
978 init_em(fplog, CG, cr, inputrec,
979 state_global, top_global, s_min, &top, &f, &f_global,
980 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
981 nfile, fnm, &outf, &mdebin);
983 /* Print to log file */
984 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
986 /* Max number of steps */
987 number_steps = inputrec->nsteps;
991 sp_header(stderr, CG, inputrec->em_tol, number_steps);
995 sp_header(fplog, CG, inputrec->em_tol, number_steps);
998 /* Call the force routine and some auxiliary (neighboursearching etc.) */
999 /* do_force always puts the charge groups in the box and shifts again
1000 * We do not unshift, so molecules are always whole in congrad.c
1002 evaluate_energy(fplog, cr,
1003 top_global, s_min, top,
1004 inputrec, nrnb, wcycle, gstat,
1005 vsite, constr, fcd, graph, mdatoms, fr,
1006 mu_tot, enerd, vir, pres, -1, TRUE);
1011 /* Copy stuff to the energy bin for easy printing etc. */
1012 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1013 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1014 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1016 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1017 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1018 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1022 /* Estimate/guess the initial stepsize */
1023 stepsize = inputrec->em_stepsize/s_min->fnorm;
1027 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1028 s_min->fmax, s_min->a_fmax+1);
1029 fprintf(stderr, " F-Norm = %12.5e\n",
1030 s_min->fnorm/sqrt(state_global->natoms));
1031 fprintf(stderr, "\n");
1032 /* and copy to the log file too... */
1033 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1034 s_min->fmax, s_min->a_fmax+1);
1035 fprintf(fplog, " F-Norm = %12.5e\n",
1036 s_min->fnorm/sqrt(state_global->natoms));
1037 fprintf(fplog, "\n");
1039 /* Start the loop over CG steps.
1040 * Each successful step is counted, and we continue until
1041 * we either converge or reach the max number of steps.
1044 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1047 /* start taking steps in a new direction
1048 * First time we enter the routine, beta=0, and the direction is
1049 * simply the negative gradient.
1052 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1057 for (i = 0; i < mdatoms->homenr; i++)
1059 if (mdatoms->cFREEZE)
1061 gf = mdatoms->cFREEZE[i];
1063 for (m = 0; m < DIM; m++)
1065 if (!inputrec->opts.nFreeze[gf][m])
1067 p[i][m] = sf[i][m] + beta*p[i][m];
1068 gpa -= p[i][m]*sf[i][m];
1069 /* f is negative gradient, thus the sign */
1078 /* Sum the gradient along the line across CPUs */
1081 gmx_sumd(1, &gpa, cr);
1084 /* Calculate the norm of the search vector */
1085 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1087 /* Just in case stepsize reaches zero due to numerical precision... */
1090 stepsize = inputrec->em_stepsize/pnorm;
1094 * Double check the value of the derivative in the search direction.
1095 * If it is positive it must be due to the old information in the
1096 * CG formula, so just remove that and start over with beta=0.
1097 * This corresponds to a steepest descent step.
1102 step--; /* Don't count this step since we are restarting */
1103 continue; /* Go back to the beginning of the big for-loop */
1106 /* Calculate minimum allowed stepsize, before the average (norm)
1107 * relative change in coordinate is smaller than precision
1110 for (i = 0; i < mdatoms->homenr; i++)
1112 for (m = 0; m < DIM; m++)
1114 tmp = fabs(s_min->s.x[i][m]);
1123 /* Add up from all CPUs */
1126 gmx_sumd(1, &minstep, cr);
1129 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1131 if (stepsize < minstep)
1137 /* Write coordinates if necessary */
1138 do_x = do_per_step(step, inputrec->nstxout);
1139 do_f = do_per_step(step, inputrec->nstfout);
1141 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1142 top_global, inputrec, step,
1143 s_min, state_global, f_global);
1145 /* Take a step downhill.
1146 * In theory, we should minimize the function along this direction.
1147 * That is quite possible, but it turns out to take 5-10 function evaluations
1148 * for each line. However, we dont really need to find the exact minimum -
1149 * it is much better to start a new CG step in a modified direction as soon
1150 * as we are close to it. This will save a lot of energy evaluations.
1152 * In practice, we just try to take a single step.
1153 * If it worked (i.e. lowered the energy), we increase the stepsize but
1154 * the continue straight to the next CG step without trying to find any minimum.
1155 * If it didn't work (higher energy), there must be a minimum somewhere between
1156 * the old position and the new one.
1158 * Due to the finite numerical accuracy, it turns out that it is a good idea
1159 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1160 * This leads to lower final energies in the tests I've done. / Erik
1162 s_a->epot = s_min->epot;
1164 c = a + stepsize; /* reference position along line is zero */
1166 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1168 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1169 s_min, top, mdatoms, fr, vsite, constr,
1173 /* Take a trial step (new coords in s_c) */
1174 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1175 constr, top, nrnb, wcycle, -1);
1178 /* Calculate energy for the trial step */
1179 evaluate_energy(fplog, cr,
1180 top_global, s_c, top,
1181 inputrec, nrnb, wcycle, gstat,
1182 vsite, constr, fcd, graph, mdatoms, fr,
1183 mu_tot, enerd, vir, pres, -1, FALSE);
1185 /* Calc derivative along line */
1189 for (i = 0; i < mdatoms->homenr; i++)
1191 for (m = 0; m < DIM; m++)
1193 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1196 /* Sum the gradient along the line across CPUs */
1199 gmx_sumd(1, &gpc, cr);
1202 /* This is the max amount of increase in energy we tolerate */
1203 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1205 /* Accept the step if the energy is lower, or if it is not significantly higher
1206 * and the line derivative is still negative.
1208 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1211 /* Great, we found a better energy. Increase step for next iteration
1212 * if we are still going down, decrease it otherwise
1216 stepsize *= 1.618034; /* The golden section */
1220 stepsize *= 0.618034; /* 1/golden section */
1225 /* New energy is the same or higher. We will have to do some work
1226 * to find a smaller value in the interval. Take smaller step next time!
1229 stepsize *= 0.618034;
1235 /* OK, if we didn't find a lower value we will have to locate one now - there must
1236 * be one in the interval [a=0,c].
1237 * The same thing is valid here, though: Don't spend dozens of iterations to find
1238 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1239 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1241 * I also have a safeguard for potentially really patological functions so we never
1242 * take more than 20 steps before we give up ...
1244 * If we already found a lower value we just skip this step and continue to the update.
1252 /* Select a new trial point.
1253 * If the derivatives at points a & c have different sign we interpolate to zero,
1254 * otherwise just do a bisection.
1256 if (gpa < 0 && gpc > 0)
1258 b = a + gpa*(a-c)/(gpc-gpa);
1265 /* safeguard if interpolation close to machine accuracy causes errors:
1266 * never go outside the interval
1268 if (b <= a || b >= c)
1273 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1275 /* Reload the old state */
1276 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1277 s_min, top, mdatoms, fr, vsite, constr,
1281 /* Take a trial step to this new point - new coords in s_b */
1282 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1283 constr, top, nrnb, wcycle, -1);
1286 /* Calculate energy for the trial step */
1287 evaluate_energy(fplog, cr,
1288 top_global, s_b, top,
1289 inputrec, nrnb, wcycle, gstat,
1290 vsite, constr, fcd, graph, mdatoms, fr,
1291 mu_tot, enerd, vir, pres, -1, FALSE);
1293 /* p does not change within a step, but since the domain decomposition
1294 * might change, we have to use cg_p of s_b here.
1299 for (i = 0; i < mdatoms->homenr; i++)
1301 for (m = 0; m < DIM; m++)
1303 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1306 /* Sum the gradient along the line across CPUs */
1309 gmx_sumd(1, &gpb, cr);
1314 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1315 s_a->epot, s_b->epot, s_c->epot, gpb);
1318 epot_repl = s_b->epot;
1320 /* Keep one of the intervals based on the value of the derivative at the new point */
1323 /* Replace c endpoint with b */
1324 swap_em_state(s_b, s_c);
1330 /* Replace a endpoint with b */
1331 swap_em_state(s_b, s_a);
1337 * Stop search as soon as we find a value smaller than the endpoints.
1338 * Never run more than 20 steps, no matter what.
1342 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1345 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1348 /* OK. We couldn't find a significantly lower energy.
1349 * If beta==0 this was steepest descent, and then we give up.
1350 * If not, set beta=0 and restart with steepest descent before quitting.
1360 /* Reset memory before giving up */
1366 /* Select min energy state of A & C, put the best in B.
1368 if (s_c->epot < s_a->epot)
1372 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1373 s_c->epot, s_a->epot);
1375 swap_em_state(s_b, s_c);
1383 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1384 s_a->epot, s_c->epot);
1386 swap_em_state(s_b, s_a);
1396 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1399 swap_em_state(s_b, s_c);
1404 /* new search direction */
1405 /* beta = 0 means forget all memory and restart with steepest descents. */
1406 if (nstcg && ((step % nstcg) == 0))
1412 /* s_min->fnorm cannot be zero, because then we would have converged
1416 /* Polak-Ribiere update.
1417 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1419 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1421 /* Limit beta to prevent oscillations */
1422 if (fabs(beta) > 5.0)
1428 /* update positions */
1429 swap_em_state(s_min, s_b);
1432 /* Print it if necessary */
1437 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1438 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1439 s_min->fmax, s_min->a_fmax+1);
1441 /* Store the new (lower) energies */
1442 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1443 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1444 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1446 do_log = do_per_step(step, inputrec->nstlog);
1447 do_ene = do_per_step(step, inputrec->nstenergy);
1450 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1452 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1453 do_log ? fplog : NULL, step, step, eprNORMAL,
1454 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1457 /* Stop when the maximum force lies below tolerance.
1458 * If we have reached machine precision, converged is already set to true.
1460 converged = converged || (s_min->fmax < inputrec->em_tol);
1462 } /* End of the loop */
1466 step--; /* we never took that last step in this case */
1469 if (s_min->fmax > inputrec->em_tol)
1473 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1474 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1481 /* If we printed energy and/or logfile last step (which was the last step)
1482 * we don't have to do it again, but otherwise print the final values.
1486 /* Write final value to log since we didn't do anything the last step */
1487 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1489 if (!do_ene || !do_log)
1491 /* Write final energy file entries */
1492 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1493 !do_log ? fplog : NULL, step, step, eprNORMAL,
1494 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1498 /* Print some stuff... */
1501 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1505 * For accurate normal mode calculation it is imperative that we
1506 * store the last conformation into the full precision binary trajectory.
1508 * However, we should only do it if we did NOT already write this step
1509 * above (which we did if do_x or do_f was true).
1511 do_x = !do_per_step(step, inputrec->nstxout);
1512 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1514 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1515 top_global, inputrec, step,
1516 s_min, state_global, f_global);
1518 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1522 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1523 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1524 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1525 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1527 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1530 finish_em(cr, outf, walltime_accounting, wcycle);
1532 /* To print the actual number of steps we needed somewhere */
1533 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1536 } /* That's all folks */
1539 double do_lbfgs(FILE *fplog, t_commrec *cr,
1540 int nfile, const t_filenm fnm[],
1541 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1542 int gmx_unused nstglobalcomm,
1543 gmx_vsite_t *vsite, gmx_constr_t constr,
1544 int gmx_unused stepout,
1545 t_inputrec *inputrec,
1546 gmx_mtop_t *top_global, t_fcdata *fcd,
1549 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1550 gmx_edsam_t gmx_unused ed,
1552 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1553 gmx_membed_t gmx_unused membed,
1554 real gmx_unused cpt_period, real gmx_unused max_hours,
1555 const char gmx_unused *deviceOptions,
1556 unsigned long gmx_unused Flags,
1557 gmx_walltime_accounting_t walltime_accounting)
1559 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1561 gmx_localtop_t *top;
1562 gmx_enerdata_t *enerd;
1564 gmx_global_stat_t gstat;
1567 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1568 double stepsize, gpa, gpb, gpc, tmp, minstep;
1569 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1570 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1571 real a, b, c, maxdelta, delta;
1572 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1573 real dgdx, dgdg, sq, yr, beta;
1575 gmx_bool converged, first;
1578 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1580 int start, end, number_steps;
1582 int i, k, m, n, nfmax, gf, step;
1589 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1594 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).");
1597 n = 3*state->natoms;
1598 nmaxcorr = inputrec->nbfgscorr;
1600 /* Allocate memory */
1601 /* Use pointers to real so we dont have to loop over both atoms and
1602 * dimensions all the time...
1603 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1604 * that point to the same memory.
1617 snew(rho, nmaxcorr);
1618 snew(alpha, nmaxcorr);
1621 for (i = 0; i < nmaxcorr; i++)
1627 for (i = 0; i < nmaxcorr; i++)
1636 init_em(fplog, LBFGS, cr, inputrec,
1637 state, top_global, &ems, &top, &f, &f_global,
1638 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1639 nfile, fnm, &outf, &mdebin);
1640 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1641 * so we free some memory again.
1646 xx = (real *)state->x;
1650 end = mdatoms->homenr;
1652 /* Print to log file */
1653 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1655 do_log = do_ene = do_x = do_f = TRUE;
1657 /* Max number of steps */
1658 number_steps = inputrec->nsteps;
1660 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1662 for (i = start; i < end; i++)
1664 if (mdatoms->cFREEZE)
1666 gf = mdatoms->cFREEZE[i];
1668 for (m = 0; m < DIM; m++)
1670 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1675 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1679 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1684 construct_vsites(vsite, state->x, 1, NULL,
1685 top->idef.iparams, top->idef.il,
1686 fr->ePBC, fr->bMolPBC, cr, state->box);
1689 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1690 /* do_force always puts the charge groups in the box and shifts again
1691 * We do not unshift, so molecules are always whole
1696 evaluate_energy(fplog, cr,
1697 top_global, &ems, top,
1698 inputrec, nrnb, wcycle, gstat,
1699 vsite, constr, fcd, graph, mdatoms, fr,
1700 mu_tot, enerd, vir, pres, -1, TRUE);
1705 /* Copy stuff to the energy bin for easy printing etc. */
1706 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1707 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1708 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1710 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1711 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1712 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1716 /* This is the starting energy */
1717 Epot = enerd->term[F_EPOT];
1723 /* Set the initial step.
1724 * since it will be multiplied by the non-normalized search direction
1725 * vector (force vector the first time), we scale it by the
1726 * norm of the force.
1731 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1732 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1733 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1734 fprintf(stderr, "\n");
1735 /* and copy to the log file too... */
1736 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1737 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1738 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1739 fprintf(fplog, "\n");
1743 for (i = 0; i < n; i++)
1747 dx[point][i] = ff[i]; /* Initial search direction */
1755 stepsize = 1.0/fnorm;
1758 /* Start the loop over BFGS steps.
1759 * Each successful step is counted, and we continue until
1760 * we either converge or reach the max number of steps.
1765 /* Set the gradient from the force */
1767 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1770 /* Write coordinates if necessary */
1771 do_x = do_per_step(step, inputrec->nstxout);
1772 do_f = do_per_step(step, inputrec->nstfout);
1777 mdof_flags |= MDOF_X;
1782 mdof_flags |= MDOF_F;
1785 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1786 top_global, step, (real)step, state, state, f, f);
1788 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1790 /* pointer to current direction - point=0 first time here */
1793 /* calculate line gradient */
1794 for (gpa = 0, i = 0; i < n; i++)
1799 /* Calculate minimum allowed stepsize, before the average (norm)
1800 * relative change in coordinate is smaller than precision
1802 for (minstep = 0, i = 0; i < n; i++)
1812 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1814 if (stepsize < minstep)
1820 /* Store old forces and coordinates */
1821 for (i = 0; i < n; i++)
1830 for (i = 0; i < n; i++)
1835 /* Take a step downhill.
1836 * In theory, we should minimize the function along this direction.
1837 * That is quite possible, but it turns out to take 5-10 function evaluations
1838 * for each line. However, we dont really need to find the exact minimum -
1839 * it is much better to start a new BFGS step in a modified direction as soon
1840 * as we are close to it. This will save a lot of energy evaluations.
1842 * In practice, we just try to take a single step.
1843 * If it worked (i.e. lowered the energy), we increase the stepsize but
1844 * the continue straight to the next BFGS step without trying to find any minimum.
1845 * If it didn't work (higher energy), there must be a minimum somewhere between
1846 * the old position and the new one.
1848 * Due to the finite numerical accuracy, it turns out that it is a good idea
1849 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1850 * This leads to lower final energies in the tests I've done. / Erik
1855 c = a + stepsize; /* reference position along line is zero */
1857 /* Check stepsize first. We do not allow displacements
1858 * larger than emstep.
1864 for (i = 0; i < n; i++)
1867 if (delta > maxdelta)
1872 if (maxdelta > inputrec->em_stepsize)
1877 while (maxdelta > inputrec->em_stepsize);
1879 /* Take a trial step */
1880 for (i = 0; i < n; i++)
1882 xc[i] = lastx[i] + c*s[i];
1886 /* Calculate energy for the trial step */
1887 ems.s.x = (rvec *)xc;
1889 evaluate_energy(fplog, cr,
1890 top_global, &ems, top,
1891 inputrec, nrnb, wcycle, gstat,
1892 vsite, constr, fcd, graph, mdatoms, fr,
1893 mu_tot, enerd, vir, pres, step, FALSE);
1896 /* Calc derivative along line */
1897 for (gpc = 0, i = 0; i < n; i++)
1899 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1901 /* Sum the gradient along the line across CPUs */
1904 gmx_sumd(1, &gpc, cr);
1907 /* This is the max amount of increase in energy we tolerate */
1908 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1910 /* Accept the step if the energy is lower, or if it is not significantly higher
1911 * and the line derivative is still negative.
1913 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1916 /* Great, we found a better energy. Increase step for next iteration
1917 * if we are still going down, decrease it otherwise
1921 stepsize *= 1.618034; /* The golden section */
1925 stepsize *= 0.618034; /* 1/golden section */
1930 /* New energy is the same or higher. We will have to do some work
1931 * to find a smaller value in the interval. Take smaller step next time!
1934 stepsize *= 0.618034;
1937 /* OK, if we didn't find a lower value we will have to locate one now - there must
1938 * be one in the interval [a=0,c].
1939 * The same thing is valid here, though: Don't spend dozens of iterations to find
1940 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1941 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1943 * I also have a safeguard for potentially really patological functions so we never
1944 * take more than 20 steps before we give up ...
1946 * If we already found a lower value we just skip this step and continue to the update.
1955 /* Select a new trial point.
1956 * If the derivatives at points a & c have different sign we interpolate to zero,
1957 * otherwise just do a bisection.
1960 if (gpa < 0 && gpc > 0)
1962 b = a + gpa*(a-c)/(gpc-gpa);
1969 /* safeguard if interpolation close to machine accuracy causes errors:
1970 * never go outside the interval
1972 if (b <= a || b >= c)
1977 /* Take a trial step */
1978 for (i = 0; i < n; i++)
1980 xb[i] = lastx[i] + b*s[i];
1984 /* Calculate energy for the trial step */
1985 ems.s.x = (rvec *)xb;
1987 evaluate_energy(fplog, cr,
1988 top_global, &ems, top,
1989 inputrec, nrnb, wcycle, gstat,
1990 vsite, constr, fcd, graph, mdatoms, fr,
1991 mu_tot, enerd, vir, pres, step, FALSE);
1996 for (gpb = 0, i = 0; i < n; i++)
1998 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2001 /* Sum the gradient along the line across CPUs */
2004 gmx_sumd(1, &gpb, cr);
2007 /* Keep one of the intervals based on the value of the derivative at the new point */
2010 /* Replace c endpoint with b */
2014 /* swap coord pointers b/c */
2024 /* Replace a endpoint with b */
2028 /* swap coord pointers a/b */
2038 * Stop search as soon as we find a value smaller than the endpoints,
2039 * or if the tolerance is below machine precision.
2040 * Never run more than 20 steps, no matter what.
2044 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2046 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2048 /* OK. We couldn't find a significantly lower energy.
2049 * If ncorr==0 this was steepest descent, and then we give up.
2050 * If not, reset memory to restart as steepest descent before quitting.
2062 /* Search in gradient direction */
2063 for (i = 0; i < n; i++)
2065 dx[point][i] = ff[i];
2067 /* Reset stepsize */
2068 stepsize = 1.0/fnorm;
2073 /* Select min energy state of A & C, put the best in xx/ff/Epot
2079 for (i = 0; i < n; i++)
2090 for (i = 0; i < n; i++)
2104 for (i = 0; i < n; i++)
2112 /* Update the memory information, and calculate a new
2113 * approximation of the inverse hessian
2116 /* Have new data in Epot, xx, ff */
2117 if (ncorr < nmaxcorr)
2122 for (i = 0; i < n; i++)
2124 dg[point][i] = lastf[i]-ff[i];
2125 dx[point][i] *= stepsize;
2130 for (i = 0; i < n; i++)
2132 dgdg += dg[point][i]*dg[point][i];
2133 dgdx += dg[point][i]*dx[point][i];
2138 rho[point] = 1.0/dgdx;
2141 if (point >= nmaxcorr)
2147 for (i = 0; i < n; i++)
2154 /* Recursive update. First go back over the memory points */
2155 for (k = 0; k < ncorr; k++)
2164 for (i = 0; i < n; i++)
2166 sq += dx[cp][i]*p[i];
2169 alpha[cp] = rho[cp]*sq;
2171 for (i = 0; i < n; i++)
2173 p[i] -= alpha[cp]*dg[cp][i];
2177 for (i = 0; i < n; i++)
2182 /* And then go forward again */
2183 for (k = 0; k < ncorr; k++)
2186 for (i = 0; i < n; i++)
2188 yr += p[i]*dg[cp][i];
2192 beta = alpha[cp]-beta;
2194 for (i = 0; i < n; i++)
2196 p[i] += beta*dx[cp][i];
2206 for (i = 0; i < n; i++)
2210 dx[point][i] = p[i];
2220 /* Test whether the convergence criterion is met */
2221 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2223 /* Print it if necessary */
2228 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2229 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2231 /* Store the new (lower) energies */
2232 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2233 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2234 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2235 do_log = do_per_step(step, inputrec->nstlog);
2236 do_ene = do_per_step(step, inputrec->nstenergy);
2239 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2241 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2242 do_log ? fplog : NULL, step, step, eprNORMAL,
2243 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2246 /* Stop when the maximum force lies below tolerance.
2247 * If we have reached machine precision, converged is already set to true.
2250 converged = converged || (fmax < inputrec->em_tol);
2252 } /* End of the loop */
2256 step--; /* we never took that last step in this case */
2259 if (fmax > inputrec->em_tol)
2263 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2264 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2269 /* If we printed energy and/or logfile last step (which was the last step)
2270 * we don't have to do it again, but otherwise print the final values.
2272 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2274 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2276 if (!do_ene || !do_log) /* Write final energy file entries */
2278 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2279 !do_log ? fplog : NULL, step, step, eprNORMAL,
2280 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2283 /* Print some stuff... */
2286 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2290 * For accurate normal mode calculation it is imperative that we
2291 * store the last conformation into the full precision binary trajectory.
2293 * However, we should only do it if we did NOT already write this step
2294 * above (which we did if do_x or do_f was true).
2296 do_x = !do_per_step(step, inputrec->nstxout);
2297 do_f = !do_per_step(step, inputrec->nstfout);
2298 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2299 top_global, inputrec, step,
2304 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2305 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2306 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2307 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2309 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2312 finish_em(cr, outf, walltime_accounting, wcycle);
2314 /* To print the actual number of steps we needed somewhere */
2315 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2318 } /* That's all folks */
2321 double do_steep(FILE *fplog, t_commrec *cr,
2322 int nfile, const t_filenm fnm[],
2323 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2324 int gmx_unused nstglobalcomm,
2325 gmx_vsite_t *vsite, gmx_constr_t constr,
2326 int gmx_unused stepout,
2327 t_inputrec *inputrec,
2328 gmx_mtop_t *top_global, t_fcdata *fcd,
2329 t_state *state_global,
2331 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2332 gmx_edsam_t gmx_unused ed,
2334 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2335 gmx_membed_t gmx_unused membed,
2336 real gmx_unused cpt_period, real gmx_unused max_hours,
2337 const char gmx_unused *deviceOptions,
2338 unsigned long gmx_unused Flags,
2339 gmx_walltime_accounting_t walltime_accounting)
2341 const char *SD = "Steepest Descents";
2342 em_state_t *s_min, *s_try;
2344 gmx_localtop_t *top;
2345 gmx_enerdata_t *enerd;
2347 gmx_global_stat_t gstat;
2349 real stepsize, constepsize;
2353 gmx_bool bDone, bAbort, do_x, do_f;
2358 int steps_accepted = 0;
2362 s_min = init_em_state();
2363 s_try = init_em_state();
2365 /* Init em and store the local state in s_try */
2366 init_em(fplog, SD, cr, inputrec,
2367 state_global, top_global, s_try, &top, &f, &f_global,
2368 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2369 nfile, fnm, &outf, &mdebin);
2371 /* Print to log file */
2372 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2374 /* Set variables for stepsize (in nm). This is the largest
2375 * step that we are going to make in any direction.
2377 ustep = inputrec->em_stepsize;
2380 /* Max number of steps */
2381 nsteps = inputrec->nsteps;
2385 /* Print to the screen */
2386 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2390 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2393 /**** HERE STARTS THE LOOP ****
2394 * count is the counter for the number of steps
2395 * bDone will be TRUE when the minimization has converged
2396 * bAbort will be TRUE when nsteps steps have been performed or when
2397 * the stepsize becomes smaller than is reasonable for machine precision
2402 while (!bDone && !bAbort)
2404 bAbort = (nsteps >= 0) && (count == nsteps);
2406 /* set new coordinates, except for first step */
2409 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2410 s_min, stepsize, s_min->f, s_try,
2411 constr, top, nrnb, wcycle, count);
2414 evaluate_energy(fplog, cr,
2415 top_global, s_try, top,
2416 inputrec, nrnb, wcycle, gstat,
2417 vsite, constr, fcd, graph, mdatoms, fr,
2418 mu_tot, enerd, vir, pres, count, count == 0);
2422 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2427 s_min->epot = s_try->epot + 1;
2430 /* Print it if necessary */
2435 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2436 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2437 (s_try->epot < s_min->epot) ? '\n' : '\r');
2440 if (s_try->epot < s_min->epot)
2442 /* Store the new (lower) energies */
2443 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2444 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2445 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2446 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2447 do_per_step(steps_accepted, inputrec->nstdisreout),
2448 do_per_step(steps_accepted, inputrec->nstorireout),
2449 fplog, count, count, eprNORMAL, TRUE,
2450 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2455 /* Now if the new energy is smaller than the previous...
2456 * or if this is the first step!
2457 * or if we did random steps!
2460 if ( (count == 0) || (s_try->epot < s_min->epot) )
2464 /* Test whether the convergence criterion is met... */
2465 bDone = (s_try->fmax < inputrec->em_tol);
2467 /* Copy the arrays for force, positions and energy */
2468 /* The 'Min' array always holds the coords and forces of the minimal
2470 swap_em_state(s_min, s_try);
2476 /* Write to trn, if necessary */
2477 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2478 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2479 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2480 top_global, inputrec, count,
2481 s_min, state_global, f_global);
2485 /* If energy is not smaller make the step smaller... */
2488 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2490 /* Reload the old state */
2491 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2492 s_min, top, mdatoms, fr, vsite, constr,
2497 /* Determine new step */
2498 stepsize = ustep/s_min->fmax;
2500 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2502 if (count == nsteps || ustep < 1e-12)
2504 if (count == nsteps || ustep < 1e-6)
2509 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2510 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2516 } /* End of the loop */
2518 /* Print some shit... */
2521 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2523 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2524 top_global, inputrec, count,
2525 s_min, state_global, f_global);
2527 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2531 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2532 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2533 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2534 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2537 finish_em(cr, outf, walltime_accounting, wcycle);
2539 /* To print the actual number of steps we needed somewhere */
2540 inputrec->nsteps = count;
2542 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2545 } /* That's all folks */
2548 double do_nm(FILE *fplog, t_commrec *cr,
2549 int nfile, const t_filenm fnm[],
2550 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2551 int gmx_unused nstglobalcomm,
2552 gmx_vsite_t *vsite, gmx_constr_t constr,
2553 int gmx_unused stepout,
2554 t_inputrec *inputrec,
2555 gmx_mtop_t *top_global, t_fcdata *fcd,
2556 t_state *state_global,
2558 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2559 gmx_edsam_t gmx_unused ed,
2561 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2562 gmx_membed_t gmx_unused membed,
2563 real gmx_unused cpt_period, real gmx_unused max_hours,
2564 const char gmx_unused *deviceOptions,
2565 unsigned long gmx_unused Flags,
2566 gmx_walltime_accounting_t walltime_accounting)
2568 const char *NM = "Normal Mode Analysis";
2570 int natoms, atom, d;
2573 gmx_localtop_t *top;
2574 gmx_enerdata_t *enerd;
2576 gmx_global_stat_t gstat;
2578 real t, t0, lambda, lam0;
2583 gmx_bool bSparse; /* use sparse matrix storage format */
2585 gmx_sparsematrix_t * sparse_matrix = NULL;
2586 real * full_matrix = NULL;
2587 em_state_t * state_work;
2589 /* added with respect to mdrun */
2590 int i, j, k, row, col;
2591 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2597 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2600 state_work = init_em_state();
2602 /* Init em and store the local state in state_minimum */
2603 init_em(fplog, NM, cr, inputrec,
2604 state_global, top_global, state_work, &top,
2606 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2607 nfile, fnm, &outf, NULL);
2609 natoms = top_global->natoms;
2617 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2618 " which MIGHT not be accurate enough for normal mode analysis.\n"
2619 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2620 " are fairly modest even if you recompile in double precision.\n\n");
2624 /* Check if we can/should use sparse storage format.
2626 * Sparse format is only useful when the Hessian itself is sparse, which it
2627 * will be when we use a cutoff.
2628 * For small systems (n<1000) it is easier to always use full matrix format, though.
2630 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2632 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2635 else if (top_global->natoms < 1000)
2637 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2642 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2648 sz = DIM*top_global->natoms;
2650 fprintf(stderr, "Allocating Hessian memory...\n\n");
2654 sparse_matrix = gmx_sparsematrix_init(sz);
2655 sparse_matrix->compressed_symmetric = TRUE;
2659 snew(full_matrix, sz*sz);
2663 /* Initial values */
2664 t0 = inputrec->init_t;
2665 lam0 = inputrec->fepvals->init_lambda;
2673 /* Write start time and temperature */
2674 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2676 /* fudge nr of steps to nr of atoms */
2677 inputrec->nsteps = natoms*2;
2681 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2682 *(top_global->name), (int)inputrec->nsteps);
2685 nnodes = cr->nnodes;
2687 /* Make evaluate_energy do a single node force calculation */
2689 evaluate_energy(fplog, cr,
2690 top_global, state_work, top,
2691 inputrec, nrnb, wcycle, gstat,
2692 vsite, constr, fcd, graph, mdatoms, fr,
2693 mu_tot, enerd, vir, pres, -1, TRUE);
2694 cr->nnodes = nnodes;
2696 /* if forces are not small, warn user */
2697 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2699 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2700 if (state_work->fmax > 1.0e-3)
2702 md_print_info(cr, fplog,
2703 "The force is probably not small enough to "
2704 "ensure that you are at a minimum.\n"
2705 "Be aware that negative eigenvalues may occur\n"
2706 "when the resulting matrix is diagonalized.\n\n");
2709 /***********************************************************
2711 * Loop over all pairs in matrix
2713 * do_force called twice. Once with positive and
2714 * once with negative displacement
2716 ************************************************************/
2718 /* Steps are divided one by one over the nodes */
2719 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2722 for (d = 0; d < DIM; d++)
2724 x_min = state_work->s.x[atom][d];
2726 state_work->s.x[atom][d] = x_min - der_range;
2728 /* Make evaluate_energy do a single node force calculation */
2730 evaluate_energy(fplog, cr,
2731 top_global, state_work, top,
2732 inputrec, nrnb, wcycle, gstat,
2733 vsite, constr, fcd, graph, mdatoms, fr,
2734 mu_tot, enerd, vir, pres, atom*2, FALSE);
2736 for (i = 0; i < natoms; i++)
2738 copy_rvec(state_work->f[i], fneg[i]);
2741 state_work->s.x[atom][d] = x_min + der_range;
2743 evaluate_energy(fplog, cr,
2744 top_global, state_work, top,
2745 inputrec, nrnb, wcycle, gstat,
2746 vsite, constr, fcd, graph, mdatoms, fr,
2747 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2748 cr->nnodes = nnodes;
2750 /* x is restored to original */
2751 state_work->s.x[atom][d] = x_min;
2753 for (j = 0; j < natoms; j++)
2755 for (k = 0; (k < DIM); k++)
2758 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2766 #define mpi_type MPI_DOUBLE
2768 #define mpi_type MPI_FLOAT
2770 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2771 cr->mpi_comm_mygroup);
2776 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2782 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2783 cr->mpi_comm_mygroup, &stat);
2788 row = (atom + node)*DIM + d;
2790 for (j = 0; j < natoms; j++)
2792 for (k = 0; k < DIM; k++)
2798 if (col >= row && dfdx[j][k] != 0.0)
2800 gmx_sparsematrix_increment_value(sparse_matrix,
2801 row, col, dfdx[j][k]);
2806 full_matrix[row*sz+col] = dfdx[j][k];
2813 if (bVerbose && fplog)
2818 /* write progress */
2819 if (MASTER(cr) && bVerbose)
2821 fprintf(stderr, "\rFinished step %d out of %d",
2822 min(atom+nnodes, natoms), natoms);
2829 fprintf(stderr, "\n\nWriting Hessian...\n");
2830 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2833 finish_em(cr, outf, walltime_accounting, wcycle);
2835 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);