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
54 #include "gmx_fatal.h"
66 #include "md_support.h"
72 #include "gmx_wallcycle.h"
73 #include "mtop_util.h"
77 #include "gmx_omp_nthreads.h"
78 #include "md_logging.h"
81 #include "gromacs/linearalgebra/mtxio.h"
82 #include "gromacs/linearalgebra/sparsematrix.h"
83 #include "gromacs/timing/walltime_accounting.h"
94 static em_state_t *init_em_state()
100 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
101 snew(ems->s.lambda, efptNR);
106 static void print_em_start(FILE *fplog,
108 gmx_walltime_accounting_t walltime_accounting,
109 gmx_wallcycle_t wcycle,
114 walltime_accounting_start(walltime_accounting);
116 sprintf(buf, "Started %s", name);
117 print_date_and_time(fplog, cr->nodeid, buf, NULL);
119 wallcycle_start(wcycle, ewcRUN);
121 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
122 gmx_wallcycle_t wcycle)
124 wallcycle_stop(wcycle, ewcRUN);
126 walltime_accounting_end(walltime_accounting);
129 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
132 fprintf(out, "%s:\n", minimizer);
133 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
134 fprintf(out, " Number of steps = %12d\n", nsteps);
137 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
143 "\nEnergy minimization reached the maximum number "
144 "of steps before the forces reached the requested "
145 "precision Fmax < %g.\n", ftol);
150 "\nEnergy minimization has stopped, but the forces have "
151 "not converged to the requested precision Fmax < %g (which "
152 "may not be possible for your system). It stopped "
153 "because the algorithm tried to make a new step whose size "
154 "was too small, or there was no change in the energy since "
155 "last step. Either way, we regard the minimization as "
156 "converged to within the available machine precision, "
157 "given your starting configuration and EM parameters.\n%s%s",
159 sizeof(real) < sizeof(double) ?
160 "\nDouble precision normally gives you higher accuracy, but "
161 "this is often not needed for preparing to run molecular "
165 "You might need to increase your constraint accuracy, or turn\n"
166 "off constraints altogether (set constraints = none in mdp file)\n" :
169 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
174 static void print_converged(FILE *fp, const char *alg, real ftol,
175 gmx_large_int_t count, gmx_bool bDone, gmx_large_int_t nsteps,
176 real epot, real fmax, int nfmax, real fnorm)
178 char buf[STEPSTRSIZE];
182 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
183 alg, ftol, gmx_step_str(count, buf));
185 else if (count < nsteps)
187 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
188 "but did not reach the requested Fmax < %g.\n",
189 alg, gmx_step_str(count, buf), ftol);
193 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
194 alg, ftol, gmx_step_str(count, buf));
198 fprintf(fp, "Potential Energy = %21.14e\n", epot);
199 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
200 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
202 fprintf(fp, "Potential Energy = %14.7e\n", epot);
203 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
204 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
208 static void get_f_norm_max(t_commrec *cr,
209 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
210 real *fnorm, real *fmax, int *a_fmax)
213 real fmax2, fmax2_0, fam;
214 int la_max, a_max, start, end, i, m, gf;
216 /* This routine finds the largest force and returns it.
217 * On parallel machines the global max is taken.
223 start = mdatoms->start;
224 end = mdatoms->homenr + start;
225 if (mdatoms->cFREEZE)
227 for (i = start; i < end; i++)
229 gf = mdatoms->cFREEZE[i];
231 for (m = 0; m < DIM; m++)
233 if (!opts->nFreeze[gf][m])
248 for (i = start; i < end; i++)
260 if (la_max >= 0 && DOMAINDECOMP(cr))
262 a_max = cr->dd->gatindex[la_max];
270 snew(sum, 2*cr->nnodes+1);
271 sum[2*cr->nodeid] = fmax2;
272 sum[2*cr->nodeid+1] = a_max;
273 sum[2*cr->nnodes] = fnorm2;
274 gmx_sumd(2*cr->nnodes+1, sum, cr);
275 fnorm2 = sum[2*cr->nnodes];
276 /* Determine the global maximum */
277 for (i = 0; i < cr->nnodes; i++)
279 if (sum[2*i] > fmax2)
282 a_max = (int)(sum[2*i+1] + 0.5);
290 *fnorm = sqrt(fnorm2);
302 static void get_state_f_norm_max(t_commrec *cr,
303 t_grpopts *opts, t_mdatoms *mdatoms,
306 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
309 void init_em(FILE *fplog, const char *title,
310 t_commrec *cr, t_inputrec *ir,
311 t_state *state_global, gmx_mtop_t *top_global,
312 em_state_t *ems, gmx_localtop_t **top,
313 rvec **f, rvec **f_global,
314 t_nrnb *nrnb, rvec mu_tot,
315 t_forcerec *fr, gmx_enerdata_t **enerd,
316 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
317 gmx_vsite_t *vsite, gmx_constr_t constr,
318 int nfile, const t_filenm fnm[],
319 gmx_mdoutf_t **outf, t_mdebin **mdebin)
321 int start, homenr, i;
326 fprintf(fplog, "Initiating %s\n", title);
329 state_global->ngtc = 0;
331 /* Initialize lambda variables */
332 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
336 if (DOMAINDECOMP(cr))
338 *top = dd_init_local_top(top_global);
340 dd_init_local_state(cr->dd, state_global, &ems->s);
344 /* Distribute the charge groups over the nodes from the master node */
345 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
346 state_global, top_global, ir,
347 &ems->s, &ems->f, mdatoms, *top,
348 fr, vsite, NULL, constr,
350 dd_store_state(cr->dd, &ems->s);
354 snew(*f_global, top_global->natoms);
364 snew(*f, top_global->natoms);
366 /* Just copy the state */
367 ems->s = *state_global;
368 snew(ems->s.x, ems->s.nalloc);
369 snew(ems->f, ems->s.nalloc);
370 for (i = 0; i < state_global->natoms; i++)
372 copy_rvec(state_global->x[i], ems->s.x[i]);
374 copy_mat(state_global->box, ems->s.box);
376 if (PAR(cr) && ir->eI != eiNM)
378 /* Initialize the particle decomposition and split the topology */
379 *top = split_system(fplog, top_global, ir, cr);
381 pd_cg_range(cr, &fr->cg0, &fr->hcg);
385 *top = gmx_mtop_generate_local_top(top_global, ir);
389 forcerec_set_excl_load(fr, *top, cr);
391 setup_bonded_threading(fr, &(*top)->idef);
393 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
395 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
404 pd_at_range(cr, &start, &homenr);
410 homenr = top_global->natoms;
412 atoms2md(top_global, ir, 0, NULL, start, homenr, mdatoms);
413 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
417 set_vsite_top(vsite, *top, mdatoms, cr);
423 if (ir->eConstrAlg == econtSHAKE &&
424 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
426 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
427 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
430 if (!DOMAINDECOMP(cr))
432 set_constraints(constr, *top, ir, mdatoms, cr);
435 if (!ir->bContinuation)
437 /* Constrain the starting coordinates */
439 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
440 ir, NULL, cr, -1, 0, mdatoms,
441 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
442 ems->s.lambda[efptFEP], &dvdl_constr,
443 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
449 *gstat = global_stat_init(ir);
452 *outf = init_mdoutf(nfile, fnm, 0, cr, ir, NULL);
455 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
460 /* Init bin for energy stuff */
461 *mdebin = init_mdebin((*outf)->fp_ene, top_global, ir, NULL);
465 calc_shifts(ems->s.box, fr->shift_vec);
468 static void finish_em(t_commrec *cr, gmx_mdoutf_t *outf,
469 gmx_walltime_accounting_t walltime_accounting,
470 gmx_wallcycle_t wcycle)
472 if (!(cr->duty & DUTY_PME))
474 /* Tell the PME only node to finish */
475 gmx_pme_send_finish(cr);
480 em_time_end(walltime_accounting, wcycle);
483 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
492 static void copy_em_coords(em_state_t *ems, t_state *state)
496 for (i = 0; (i < state->natoms); i++)
498 copy_rvec(ems->s.x[i], state->x[i]);
502 static void write_em_traj(FILE *fplog, t_commrec *cr,
504 gmx_bool bX, gmx_bool bF, const char *confout,
505 gmx_mtop_t *top_global,
506 t_inputrec *ir, gmx_large_int_t step,
508 t_state *state_global, rvec *f_global)
512 if ((bX || bF || confout != NULL) && !DOMAINDECOMP(cr))
514 copy_em_coords(state, state_global);
521 mdof_flags |= MDOF_X;
525 mdof_flags |= MDOF_F;
527 write_traj(fplog, cr, outf, mdof_flags,
528 top_global, step, (double)step,
529 &state->s, state_global, state->f, f_global, NULL, NULL);
531 if (confout != NULL && MASTER(cr))
533 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
535 /* Make molecules whole only for confout writing */
536 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
540 write_sto_conf_mtop(confout,
541 *top_global->name, top_global,
542 state_global->x, NULL, ir->ePBC, state_global->box);
546 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
548 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
549 gmx_constr_t constr, gmx_localtop_t *top,
550 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
551 gmx_large_int_t count)
563 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
565 gmx_incons("state mismatch in do_em_step");
568 s2->flags = s1->flags;
570 if (s2->nalloc != s1->nalloc)
572 s2->nalloc = s1->nalloc;
573 srenew(s2->x, s1->nalloc);
574 srenew(ems2->f, s1->nalloc);
575 if (s2->flags & (1<<estCGP))
577 srenew(s2->cg_p, s1->nalloc);
581 s2->natoms = s1->natoms;
582 copy_mat(s1->box, s2->box);
583 /* Copy free energy state */
584 for (i = 0; i < efptNR; i++)
586 s2->lambda[i] = s1->lambda[i];
588 copy_mat(s1->box, s2->box);
591 end = md->start + md->homenr;
596 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
601 #pragma omp for schedule(static) nowait
602 for (i = start; i < end; i++)
608 for (m = 0; m < DIM; m++)
610 if (ir->opts.nFreeze[gf][m])
616 x2[i][m] = x1[i][m] + a*f[i][m];
621 if (s2->flags & (1<<estCGP))
623 /* Copy the CG p vector */
626 #pragma omp for schedule(static) nowait
627 for (i = start; i < end; i++)
629 copy_rvec(x1[i], x2[i]);
633 if (DOMAINDECOMP(cr))
635 s2->ddp_count = s1->ddp_count;
636 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
639 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
640 srenew(s2->cg_gl, s2->cg_gl_nalloc);
643 s2->ncg_gl = s1->ncg_gl;
644 #pragma omp for schedule(static) nowait
645 for (i = 0; i < s2->ncg_gl; i++)
647 s2->cg_gl[i] = s1->cg_gl[i];
649 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
655 wallcycle_start(wcycle, ewcCONSTR);
657 constrain(NULL, TRUE, TRUE, constr, &top->idef,
658 ir, NULL, cr, count, 0, md,
659 s1->x, s2->x, NULL, bMolPBC, s2->box,
660 s2->lambda[efptBONDED], &dvdl_constr,
661 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
662 wallcycle_stop(wcycle, ewcCONSTR);
666 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
667 gmx_mtop_t *top_global, t_inputrec *ir,
668 em_state_t *ems, gmx_localtop_t *top,
669 t_mdatoms *mdatoms, t_forcerec *fr,
670 gmx_vsite_t *vsite, gmx_constr_t constr,
671 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
673 /* Repartition the domain decomposition */
674 wallcycle_start(wcycle, ewcDOMDEC);
675 dd_partition_system(fplog, step, cr, FALSE, 1,
676 NULL, top_global, ir,
678 mdatoms, top, fr, vsite, NULL, constr,
679 nrnb, wcycle, FALSE);
680 dd_store_state(cr->dd, &ems->s);
681 wallcycle_stop(wcycle, ewcDOMDEC);
684 static void evaluate_energy(FILE *fplog, t_commrec *cr,
685 gmx_mtop_t *top_global,
686 em_state_t *ems, gmx_localtop_t *top,
687 t_inputrec *inputrec,
688 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
689 gmx_global_stat_t gstat,
690 gmx_vsite_t *vsite, gmx_constr_t constr,
692 t_graph *graph, t_mdatoms *mdatoms,
693 t_forcerec *fr, rvec mu_tot,
694 gmx_enerdata_t *enerd, tensor vir, tensor pres,
695 gmx_large_int_t count, gmx_bool bFirst)
700 tensor force_vir, shake_vir, ekin;
701 real dvdl_constr, prescorr, enercorr, dvdlcorr;
704 /* Set the time to the initial time, the time does not change during EM */
705 t = inputrec->init_t;
708 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
710 /* This the first state or an old state used before the last ns */
716 if (inputrec->nstlist > 0)
720 else if (inputrec->nstlist == -1)
722 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
725 gmx_sumi(1, &nabnsb, cr);
733 construct_vsites(vsite, ems->s.x, 1, NULL,
734 top->idef.iparams, top->idef.il,
735 fr->ePBC, fr->bMolPBC, graph, cr, ems->s.box);
738 if (DOMAINDECOMP(cr))
742 /* Repartition the domain decomposition */
743 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
744 ems, top, mdatoms, fr, vsite, constr,
749 /* Calc force & energy on new trial position */
750 /* do_force always puts the charge groups in the box and shifts again
751 * We do not unshift, so molecules are always whole in congrad.c
753 do_force(fplog, cr, inputrec,
754 count, nrnb, wcycle, top, &top_global->groups,
755 ems->s.box, ems->s.x, &ems->s.hist,
756 ems->f, force_vir, mdatoms, enerd, fcd,
757 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
758 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
759 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
760 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
762 /* Clear the unused shake virial and pressure */
763 clear_mat(shake_vir);
766 /* Communicate stuff when parallel */
767 if (PAR(cr) && inputrec->eI != eiNM)
769 wallcycle_start(wcycle, ewcMoveE);
771 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
772 inputrec, NULL, NULL, NULL, 1, &terminate,
773 top_global, &ems->s, FALSE,
779 wallcycle_stop(wcycle, ewcMoveE);
782 /* Calculate long range corrections to pressure and energy */
783 calc_dispcorr(fplog, inputrec, fr, count, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
784 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
785 enerd->term[F_DISPCORR] = enercorr;
786 enerd->term[F_EPOT] += enercorr;
787 enerd->term[F_PRES] += prescorr;
788 enerd->term[F_DVDL] += dvdlcorr;
790 ems->epot = enerd->term[F_EPOT];
794 /* Project out the constraint components of the force */
795 wallcycle_start(wcycle, ewcCONSTR);
797 constrain(NULL, FALSE, FALSE, constr, &top->idef,
798 inputrec, NULL, cr, count, 0, mdatoms,
799 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
800 ems->s.lambda[efptBONDED], &dvdl_constr,
801 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
802 if (fr->bSepDVDL && fplog)
804 gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
806 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
807 m_add(force_vir, shake_vir, vir);
808 wallcycle_stop(wcycle, ewcCONSTR);
812 copy_mat(force_vir, vir);
816 enerd->term[F_PRES] =
817 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
819 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
821 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
823 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
827 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
829 em_state_t *s_min, em_state_t *s_b)
833 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
835 unsigned char *grpnrFREEZE;
839 fprintf(debug, "Doing reorder_partsum\n");
845 cgs_gl = dd_charge_groups_global(cr->dd);
846 index = cgs_gl->index;
848 /* Collect fm in a global vector fmg.
849 * This conflicts with the spirit of domain decomposition,
850 * but to fully optimize this a much more complicated algorithm is required.
852 snew(fmg, mtop->natoms);
854 ncg = s_min->s.ncg_gl;
855 cg_gl = s_min->s.cg_gl;
857 for (c = 0; c < ncg; c++)
862 for (a = a0; a < a1; a++)
864 copy_rvec(fm[i], fmg[a]);
868 gmx_sum(mtop->natoms*3, fmg[0], cr);
870 /* Now we will determine the part of the sum for the cgs in state s_b */
872 cg_gl = s_b->s.cg_gl;
876 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
877 for (c = 0; c < ncg; c++)
882 for (a = a0; a < a1; a++)
884 if (mdatoms->cFREEZE && grpnrFREEZE)
888 for (m = 0; m < DIM; m++)
890 if (!opts->nFreeze[gf][m])
892 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
904 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
906 em_state_t *s_min, em_state_t *s_b)
912 /* This is just the classical Polak-Ribiere calculation of beta;
913 * it looks a bit complicated since we take freeze groups into account,
914 * and might have to sum it in parallel runs.
917 if (!DOMAINDECOMP(cr) ||
918 (s_min->s.ddp_count == cr->dd->ddp_count &&
919 s_b->s.ddp_count == cr->dd->ddp_count))
925 /* This part of code can be incorrect with DD,
926 * since the atom ordering in s_b and s_min might differ.
928 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
930 if (mdatoms->cFREEZE)
932 gf = mdatoms->cFREEZE[i];
934 for (m = 0; m < DIM; m++)
936 if (!opts->nFreeze[gf][m])
938 sum += (fb[i][m] - fm[i][m])*fb[i][m];
945 /* We need to reorder cgs while summing */
946 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
950 gmx_sumd(1, &sum, cr);
953 return sum/sqr(s_min->fnorm);
956 double do_cg(FILE *fplog, t_commrec *cr,
957 int nfile, const t_filenm fnm[],
958 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
959 int gmx_unused nstglobalcomm,
960 gmx_vsite_t *vsite, gmx_constr_t constr,
961 int gmx_unused stepout,
962 t_inputrec *inputrec,
963 gmx_mtop_t *top_global, t_fcdata *fcd,
964 t_state *state_global,
966 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
967 gmx_edsam_t gmx_unused ed,
969 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
970 gmx_membed_t gmx_unused membed,
971 real gmx_unused cpt_period, real gmx_unused max_hours,
972 const char gmx_unused *deviceOptions,
973 unsigned long gmx_unused Flags,
974 gmx_walltime_accounting_t walltime_accounting)
976 const char *CG = "Polak-Ribiere Conjugate Gradients";
978 em_state_t *s_min, *s_a, *s_b, *s_c;
980 gmx_enerdata_t *enerd;
982 gmx_global_stat_t gstat;
984 rvec *f_global, *p, *sf, *sfm;
985 double gpa, gpb, gpc, tmp, sum[2], minstep;
988 real a, b, c, beta = 0.0;
992 gmx_bool converged, foundlower;
994 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
996 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
998 int i, m, gf, step, nminstep;
1003 s_min = init_em_state();
1004 s_a = init_em_state();
1005 s_b = init_em_state();
1006 s_c = init_em_state();
1008 /* Init em and store the local state in s_min */
1009 init_em(fplog, CG, cr, inputrec,
1010 state_global, top_global, s_min, &top, &f, &f_global,
1011 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1012 nfile, fnm, &outf, &mdebin);
1014 /* Print to log file */
1015 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1017 /* Max number of steps */
1018 number_steps = inputrec->nsteps;
1022 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1026 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1029 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1030 /* do_force always puts the charge groups in the box and shifts again
1031 * We do not unshift, so molecules are always whole in congrad.c
1033 evaluate_energy(fplog, cr,
1034 top_global, s_min, top,
1035 inputrec, nrnb, wcycle, gstat,
1036 vsite, constr, fcd, graph, mdatoms, fr,
1037 mu_tot, enerd, vir, pres, -1, TRUE);
1042 /* Copy stuff to the energy bin for easy printing etc. */
1043 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1044 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1045 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1047 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1048 print_ebin(outf->fp_ene, TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1049 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1053 /* Estimate/guess the initial stepsize */
1054 stepsize = inputrec->em_stepsize/s_min->fnorm;
1058 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1059 s_min->fmax, s_min->a_fmax+1);
1060 fprintf(stderr, " F-Norm = %12.5e\n",
1061 s_min->fnorm/sqrt(state_global->natoms));
1062 fprintf(stderr, "\n");
1063 /* and copy to the log file too... */
1064 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1065 s_min->fmax, s_min->a_fmax+1);
1066 fprintf(fplog, " F-Norm = %12.5e\n",
1067 s_min->fnorm/sqrt(state_global->natoms));
1068 fprintf(fplog, "\n");
1070 /* Start the loop over CG steps.
1071 * Each successful step is counted, and we continue until
1072 * we either converge or reach the max number of steps.
1075 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1078 /* start taking steps in a new direction
1079 * First time we enter the routine, beta=0, and the direction is
1080 * simply the negative gradient.
1083 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1088 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1090 if (mdatoms->cFREEZE)
1092 gf = mdatoms->cFREEZE[i];
1094 for (m = 0; m < DIM; m++)
1096 if (!inputrec->opts.nFreeze[gf][m])
1098 p[i][m] = sf[i][m] + beta*p[i][m];
1099 gpa -= p[i][m]*sf[i][m];
1100 /* f is negative gradient, thus the sign */
1109 /* Sum the gradient along the line across CPUs */
1112 gmx_sumd(1, &gpa, cr);
1115 /* Calculate the norm of the search vector */
1116 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1118 /* Just in case stepsize reaches zero due to numerical precision... */
1121 stepsize = inputrec->em_stepsize/pnorm;
1125 * Double check the value of the derivative in the search direction.
1126 * If it is positive it must be due to the old information in the
1127 * CG formula, so just remove that and start over with beta=0.
1128 * This corresponds to a steepest descent step.
1133 step--; /* Don't count this step since we are restarting */
1134 continue; /* Go back to the beginning of the big for-loop */
1137 /* Calculate minimum allowed stepsize, before the average (norm)
1138 * relative change in coordinate is smaller than precision
1141 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1143 for (m = 0; m < DIM; m++)
1145 tmp = fabs(s_min->s.x[i][m]);
1154 /* Add up from all CPUs */
1157 gmx_sumd(1, &minstep, cr);
1160 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1162 if (stepsize < minstep)
1168 /* Write coordinates if necessary */
1169 do_x = do_per_step(step, inputrec->nstxout);
1170 do_f = do_per_step(step, inputrec->nstfout);
1172 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1173 top_global, inputrec, step,
1174 s_min, state_global, f_global);
1176 /* Take a step downhill.
1177 * In theory, we should minimize the function along this direction.
1178 * That is quite possible, but it turns out to take 5-10 function evaluations
1179 * for each line. However, we dont really need to find the exact minimum -
1180 * it is much better to start a new CG step in a modified direction as soon
1181 * as we are close to it. This will save a lot of energy evaluations.
1183 * In practice, we just try to take a single step.
1184 * If it worked (i.e. lowered the energy), we increase the stepsize but
1185 * the continue straight to the next CG step without trying to find any minimum.
1186 * If it didn't work (higher energy), there must be a minimum somewhere between
1187 * the old position and the new one.
1189 * Due to the finite numerical accuracy, it turns out that it is a good idea
1190 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1191 * This leads to lower final energies in the tests I've done. / Erik
1193 s_a->epot = s_min->epot;
1195 c = a + stepsize; /* reference position along line is zero */
1197 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1199 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1200 s_min, top, mdatoms, fr, vsite, constr,
1204 /* Take a trial step (new coords in s_c) */
1205 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1206 constr, top, nrnb, wcycle, -1);
1209 /* Calculate energy for the trial step */
1210 evaluate_energy(fplog, cr,
1211 top_global, s_c, top,
1212 inputrec, nrnb, wcycle, gstat,
1213 vsite, constr, fcd, graph, mdatoms, fr,
1214 mu_tot, enerd, vir, pres, -1, FALSE);
1216 /* Calc derivative along line */
1220 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1222 for (m = 0; m < DIM; m++)
1224 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1227 /* Sum the gradient along the line across CPUs */
1230 gmx_sumd(1, &gpc, cr);
1233 /* This is the max amount of increase in energy we tolerate */
1234 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1236 /* Accept the step if the energy is lower, or if it is not significantly higher
1237 * and the line derivative is still negative.
1239 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1242 /* Great, we found a better energy. Increase step for next iteration
1243 * if we are still going down, decrease it otherwise
1247 stepsize *= 1.618034; /* The golden section */
1251 stepsize *= 0.618034; /* 1/golden section */
1256 /* New energy is the same or higher. We will have to do some work
1257 * to find a smaller value in the interval. Take smaller step next time!
1260 stepsize *= 0.618034;
1266 /* OK, if we didn't find a lower value we will have to locate one now - there must
1267 * be one in the interval [a=0,c].
1268 * The same thing is valid here, though: Don't spend dozens of iterations to find
1269 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1270 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1272 * I also have a safeguard for potentially really patological functions so we never
1273 * take more than 20 steps before we give up ...
1275 * If we already found a lower value we just skip this step and continue to the update.
1283 /* Select a new trial point.
1284 * If the derivatives at points a & c have different sign we interpolate to zero,
1285 * otherwise just do a bisection.
1287 if (gpa < 0 && gpc > 0)
1289 b = a + gpa*(a-c)/(gpc-gpa);
1296 /* safeguard if interpolation close to machine accuracy causes errors:
1297 * never go outside the interval
1299 if (b <= a || b >= c)
1304 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1306 /* Reload the old state */
1307 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1308 s_min, top, mdatoms, fr, vsite, constr,
1312 /* Take a trial step to this new point - new coords in s_b */
1313 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1314 constr, top, nrnb, wcycle, -1);
1317 /* Calculate energy for the trial step */
1318 evaluate_energy(fplog, cr,
1319 top_global, s_b, top,
1320 inputrec, nrnb, wcycle, gstat,
1321 vsite, constr, fcd, graph, mdatoms, fr,
1322 mu_tot, enerd, vir, pres, -1, FALSE);
1324 /* p does not change within a step, but since the domain decomposition
1325 * might change, we have to use cg_p of s_b here.
1330 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1332 for (m = 0; m < DIM; m++)
1334 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1337 /* Sum the gradient along the line across CPUs */
1340 gmx_sumd(1, &gpb, cr);
1345 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1346 s_a->epot, s_b->epot, s_c->epot, gpb);
1349 epot_repl = s_b->epot;
1351 /* Keep one of the intervals based on the value of the derivative at the new point */
1354 /* Replace c endpoint with b */
1355 swap_em_state(s_b, s_c);
1361 /* Replace a endpoint with b */
1362 swap_em_state(s_b, s_a);
1368 * Stop search as soon as we find a value smaller than the endpoints.
1369 * Never run more than 20 steps, no matter what.
1373 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1376 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1379 /* OK. We couldn't find a significantly lower energy.
1380 * If beta==0 this was steepest descent, and then we give up.
1381 * If not, set beta=0 and restart with steepest descent before quitting.
1391 /* Reset memory before giving up */
1397 /* Select min energy state of A & C, put the best in B.
1399 if (s_c->epot < s_a->epot)
1403 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1404 s_c->epot, s_a->epot);
1406 swap_em_state(s_b, s_c);
1414 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1415 s_a->epot, s_c->epot);
1417 swap_em_state(s_b, s_a);
1427 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1430 swap_em_state(s_b, s_c);
1435 /* new search direction */
1436 /* beta = 0 means forget all memory and restart with steepest descents. */
1437 if (nstcg && ((step % nstcg) == 0))
1443 /* s_min->fnorm cannot be zero, because then we would have converged
1447 /* Polak-Ribiere update.
1448 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1450 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1452 /* Limit beta to prevent oscillations */
1453 if (fabs(beta) > 5.0)
1459 /* update positions */
1460 swap_em_state(s_min, s_b);
1463 /* Print it if necessary */
1468 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1469 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1470 s_min->fmax, s_min->a_fmax+1);
1472 /* Store the new (lower) energies */
1473 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1474 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1475 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1477 do_log = do_per_step(step, inputrec->nstlog);
1478 do_ene = do_per_step(step, inputrec->nstenergy);
1481 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1483 print_ebin(outf->fp_ene, do_ene, FALSE, FALSE,
1484 do_log ? fplog : NULL, step, step, eprNORMAL,
1485 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1488 /* Stop when the maximum force lies below tolerance.
1489 * If we have reached machine precision, converged is already set to true.
1491 converged = converged || (s_min->fmax < inputrec->em_tol);
1493 } /* End of the loop */
1497 step--; /* we never took that last step in this case */
1500 if (s_min->fmax > inputrec->em_tol)
1504 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1505 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1512 /* If we printed energy and/or logfile last step (which was the last step)
1513 * we don't have to do it again, but otherwise print the final values.
1517 /* Write final value to log since we didn't do anything the last step */
1518 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1520 if (!do_ene || !do_log)
1522 /* Write final energy file entries */
1523 print_ebin(outf->fp_ene, !do_ene, FALSE, FALSE,
1524 !do_log ? fplog : NULL, step, step, eprNORMAL,
1525 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1529 /* Print some stuff... */
1532 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1536 * For accurate normal mode calculation it is imperative that we
1537 * store the last conformation into the full precision binary trajectory.
1539 * However, we should only do it if we did NOT already write this step
1540 * above (which we did if do_x or do_f was true).
1542 do_x = !do_per_step(step, inputrec->nstxout);
1543 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1545 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1546 top_global, inputrec, step,
1547 s_min, state_global, f_global);
1549 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1553 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1554 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1555 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1556 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1558 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1561 finish_em(cr, outf, walltime_accounting, wcycle);
1563 /* To print the actual number of steps we needed somewhere */
1564 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1567 } /* That's all folks */
1570 double do_lbfgs(FILE *fplog, t_commrec *cr,
1571 int nfile, const t_filenm fnm[],
1572 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1573 int gmx_unused nstglobalcomm,
1574 gmx_vsite_t *vsite, gmx_constr_t constr,
1575 int gmx_unused stepout,
1576 t_inputrec *inputrec,
1577 gmx_mtop_t *top_global, t_fcdata *fcd,
1580 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1581 gmx_edsam_t gmx_unused ed,
1583 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1584 gmx_membed_t gmx_unused membed,
1585 real gmx_unused cpt_period, real gmx_unused max_hours,
1586 const char gmx_unused *deviceOptions,
1587 unsigned long gmx_unused Flags,
1588 gmx_walltime_accounting_t walltime_accounting)
1590 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1592 gmx_localtop_t *top;
1593 gmx_enerdata_t *enerd;
1595 gmx_global_stat_t gstat;
1598 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1599 double stepsize, gpa, gpb, gpc, tmp, minstep;
1600 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1601 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1602 real a, b, c, maxdelta, delta;
1603 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1604 real dgdx, dgdg, sq, yr, beta;
1606 gmx_bool converged, first;
1609 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1611 int start, end, number_steps;
1613 int i, k, m, n, nfmax, gf, step;
1620 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1625 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).");
1628 n = 3*state->natoms;
1629 nmaxcorr = inputrec->nbfgscorr;
1631 /* Allocate memory */
1632 /* Use pointers to real so we dont have to loop over both atoms and
1633 * dimensions all the time...
1634 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1635 * that point to the same memory.
1648 snew(rho, nmaxcorr);
1649 snew(alpha, nmaxcorr);
1652 for (i = 0; i < nmaxcorr; i++)
1658 for (i = 0; i < nmaxcorr; i++)
1667 init_em(fplog, LBFGS, cr, inputrec,
1668 state, top_global, &ems, &top, &f, &f_global,
1669 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1670 nfile, fnm, &outf, &mdebin);
1671 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1672 * so we free some memory again.
1677 xx = (real *)state->x;
1680 start = mdatoms->start;
1681 end = mdatoms->homenr + start;
1683 /* Print to log file */
1684 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1686 do_log = do_ene = do_x = do_f = TRUE;
1688 /* Max number of steps */
1689 number_steps = inputrec->nsteps;
1691 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1693 for (i = start; i < end; i++)
1695 if (mdatoms->cFREEZE)
1697 gf = mdatoms->cFREEZE[i];
1699 for (m = 0; m < DIM; m++)
1701 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1706 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1710 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1715 construct_vsites(vsite, state->x, 1, NULL,
1716 top->idef.iparams, top->idef.il,
1717 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1720 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1721 /* do_force always puts the charge groups in the box and shifts again
1722 * We do not unshift, so molecules are always whole
1727 evaluate_energy(fplog, cr,
1728 top_global, &ems, top,
1729 inputrec, nrnb, wcycle, gstat,
1730 vsite, constr, fcd, graph, mdatoms, fr,
1731 mu_tot, enerd, vir, pres, -1, TRUE);
1736 /* Copy stuff to the energy bin for easy printing etc. */
1737 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1738 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1739 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1741 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1742 print_ebin(outf->fp_ene, TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1743 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1747 /* This is the starting energy */
1748 Epot = enerd->term[F_EPOT];
1754 /* Set the initial step.
1755 * since it will be multiplied by the non-normalized search direction
1756 * vector (force vector the first time), we scale it by the
1757 * norm of the force.
1762 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1763 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1764 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1765 fprintf(stderr, "\n");
1766 /* and copy to the log file too... */
1767 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1768 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1769 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1770 fprintf(fplog, "\n");
1774 for (i = 0; i < n; i++)
1778 dx[point][i] = ff[i]; /* Initial search direction */
1786 stepsize = 1.0/fnorm;
1789 /* Start the loop over BFGS steps.
1790 * Each successful step is counted, and we continue until
1791 * we either converge or reach the max number of steps.
1796 /* Set the gradient from the force */
1798 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1801 /* Write coordinates if necessary */
1802 do_x = do_per_step(step, inputrec->nstxout);
1803 do_f = do_per_step(step, inputrec->nstfout);
1808 mdof_flags |= MDOF_X;
1813 mdof_flags |= MDOF_F;
1816 write_traj(fplog, cr, outf, mdof_flags,
1817 top_global, step, (real)step, state, state, f, f, NULL, NULL);
1819 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1821 /* pointer to current direction - point=0 first time here */
1824 /* calculate line gradient */
1825 for (gpa = 0, i = 0; i < n; i++)
1830 /* Calculate minimum allowed stepsize, before the average (norm)
1831 * relative change in coordinate is smaller than precision
1833 for (minstep = 0, i = 0; i < n; i++)
1843 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1845 if (stepsize < minstep)
1851 /* Store old forces and coordinates */
1852 for (i = 0; i < n; i++)
1861 for (i = 0; i < n; i++)
1866 /* Take a step downhill.
1867 * In theory, we should minimize the function along this direction.
1868 * That is quite possible, but it turns out to take 5-10 function evaluations
1869 * for each line. However, we dont really need to find the exact minimum -
1870 * it is much better to start a new BFGS step in a modified direction as soon
1871 * as we are close to it. This will save a lot of energy evaluations.
1873 * In practice, we just try to take a single step.
1874 * If it worked (i.e. lowered the energy), we increase the stepsize but
1875 * the continue straight to the next BFGS step without trying to find any minimum.
1876 * If it didn't work (higher energy), there must be a minimum somewhere between
1877 * the old position and the new one.
1879 * Due to the finite numerical accuracy, it turns out that it is a good idea
1880 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1881 * This leads to lower final energies in the tests I've done. / Erik
1886 c = a + stepsize; /* reference position along line is zero */
1888 /* Check stepsize first. We do not allow displacements
1889 * larger than emstep.
1895 for (i = 0; i < n; i++)
1898 if (delta > maxdelta)
1903 if (maxdelta > inputrec->em_stepsize)
1908 while (maxdelta > inputrec->em_stepsize);
1910 /* Take a trial step */
1911 for (i = 0; i < n; i++)
1913 xc[i] = lastx[i] + c*s[i];
1917 /* Calculate energy for the trial step */
1918 ems.s.x = (rvec *)xc;
1920 evaluate_energy(fplog, cr,
1921 top_global, &ems, top,
1922 inputrec, nrnb, wcycle, gstat,
1923 vsite, constr, fcd, graph, mdatoms, fr,
1924 mu_tot, enerd, vir, pres, step, FALSE);
1927 /* Calc derivative along line */
1928 for (gpc = 0, i = 0; i < n; i++)
1930 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1932 /* Sum the gradient along the line across CPUs */
1935 gmx_sumd(1, &gpc, cr);
1938 /* This is the max amount of increase in energy we tolerate */
1939 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1941 /* Accept the step if the energy is lower, or if it is not significantly higher
1942 * and the line derivative is still negative.
1944 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1947 /* Great, we found a better energy. Increase step for next iteration
1948 * if we are still going down, decrease it otherwise
1952 stepsize *= 1.618034; /* The golden section */
1956 stepsize *= 0.618034; /* 1/golden section */
1961 /* New energy is the same or higher. We will have to do some work
1962 * to find a smaller value in the interval. Take smaller step next time!
1965 stepsize *= 0.618034;
1968 /* OK, if we didn't find a lower value we will have to locate one now - there must
1969 * be one in the interval [a=0,c].
1970 * The same thing is valid here, though: Don't spend dozens of iterations to find
1971 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1972 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1974 * I also have a safeguard for potentially really patological functions so we never
1975 * take more than 20 steps before we give up ...
1977 * If we already found a lower value we just skip this step and continue to the update.
1986 /* Select a new trial point.
1987 * If the derivatives at points a & c have different sign we interpolate to zero,
1988 * otherwise just do a bisection.
1991 if (gpa < 0 && gpc > 0)
1993 b = a + gpa*(a-c)/(gpc-gpa);
2000 /* safeguard if interpolation close to machine accuracy causes errors:
2001 * never go outside the interval
2003 if (b <= a || b >= c)
2008 /* Take a trial step */
2009 for (i = 0; i < n; i++)
2011 xb[i] = lastx[i] + b*s[i];
2015 /* Calculate energy for the trial step */
2016 ems.s.x = (rvec *)xb;
2018 evaluate_energy(fplog, cr,
2019 top_global, &ems, top,
2020 inputrec, nrnb, wcycle, gstat,
2021 vsite, constr, fcd, graph, mdatoms, fr,
2022 mu_tot, enerd, vir, pres, step, FALSE);
2027 for (gpb = 0, i = 0; i < n; i++)
2029 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2032 /* Sum the gradient along the line across CPUs */
2035 gmx_sumd(1, &gpb, cr);
2038 /* Keep one of the intervals based on the value of the derivative at the new point */
2041 /* Replace c endpoint with b */
2045 /* swap coord pointers b/c */
2055 /* Replace a endpoint with b */
2059 /* swap coord pointers a/b */
2069 * Stop search as soon as we find a value smaller than the endpoints,
2070 * or if the tolerance is below machine precision.
2071 * Never run more than 20 steps, no matter what.
2075 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2077 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2079 /* OK. We couldn't find a significantly lower energy.
2080 * If ncorr==0 this was steepest descent, and then we give up.
2081 * If not, reset memory to restart as steepest descent before quitting.
2093 /* Search in gradient direction */
2094 for (i = 0; i < n; i++)
2096 dx[point][i] = ff[i];
2098 /* Reset stepsize */
2099 stepsize = 1.0/fnorm;
2104 /* Select min energy state of A & C, put the best in xx/ff/Epot
2110 for (i = 0; i < n; i++)
2121 for (i = 0; i < n; i++)
2135 for (i = 0; i < n; i++)
2143 /* Update the memory information, and calculate a new
2144 * approximation of the inverse hessian
2147 /* Have new data in Epot, xx, ff */
2148 if (ncorr < nmaxcorr)
2153 for (i = 0; i < n; i++)
2155 dg[point][i] = lastf[i]-ff[i];
2156 dx[point][i] *= stepsize;
2161 for (i = 0; i < n; i++)
2163 dgdg += dg[point][i]*dg[point][i];
2164 dgdx += dg[point][i]*dx[point][i];
2169 rho[point] = 1.0/dgdx;
2172 if (point >= nmaxcorr)
2178 for (i = 0; i < n; i++)
2185 /* Recursive update. First go back over the memory points */
2186 for (k = 0; k < ncorr; k++)
2195 for (i = 0; i < n; i++)
2197 sq += dx[cp][i]*p[i];
2200 alpha[cp] = rho[cp]*sq;
2202 for (i = 0; i < n; i++)
2204 p[i] -= alpha[cp]*dg[cp][i];
2208 for (i = 0; i < n; i++)
2213 /* And then go forward again */
2214 for (k = 0; k < ncorr; k++)
2217 for (i = 0; i < n; i++)
2219 yr += p[i]*dg[cp][i];
2223 beta = alpha[cp]-beta;
2225 for (i = 0; i < n; i++)
2227 p[i] += beta*dx[cp][i];
2237 for (i = 0; i < n; i++)
2241 dx[point][i] = p[i];
2251 /* Test whether the convergence criterion is met */
2252 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2254 /* Print it if necessary */
2259 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2260 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2262 /* Store the new (lower) energies */
2263 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2264 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2265 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2266 do_log = do_per_step(step, inputrec->nstlog);
2267 do_ene = do_per_step(step, inputrec->nstenergy);
2270 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2272 print_ebin(outf->fp_ene, do_ene, FALSE, FALSE,
2273 do_log ? fplog : NULL, step, step, eprNORMAL,
2274 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2277 /* Stop when the maximum force lies below tolerance.
2278 * If we have reached machine precision, converged is already set to true.
2281 converged = converged || (fmax < inputrec->em_tol);
2283 } /* End of the loop */
2287 step--; /* we never took that last step in this case */
2290 if (fmax > inputrec->em_tol)
2294 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2295 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2300 /* If we printed energy and/or logfile last step (which was the last step)
2301 * we don't have to do it again, but otherwise print the final values.
2303 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2305 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2307 if (!do_ene || !do_log) /* Write final energy file entries */
2309 print_ebin(outf->fp_ene, !do_ene, FALSE, FALSE,
2310 !do_log ? fplog : NULL, step, step, eprNORMAL,
2311 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2314 /* Print some stuff... */
2317 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2321 * For accurate normal mode calculation it is imperative that we
2322 * store the last conformation into the full precision binary trajectory.
2324 * However, we should only do it if we did NOT already write this step
2325 * above (which we did if do_x or do_f was true).
2327 do_x = !do_per_step(step, inputrec->nstxout);
2328 do_f = !do_per_step(step, inputrec->nstfout);
2329 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2330 top_global, inputrec, step,
2335 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2336 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2337 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2338 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2340 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2343 finish_em(cr, outf, walltime_accounting, wcycle);
2345 /* To print the actual number of steps we needed somewhere */
2346 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2349 } /* That's all folks */
2352 double do_steep(FILE *fplog, t_commrec *cr,
2353 int nfile, const t_filenm fnm[],
2354 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2355 int gmx_unused nstglobalcomm,
2356 gmx_vsite_t *vsite, gmx_constr_t constr,
2357 int gmx_unused stepout,
2358 t_inputrec *inputrec,
2359 gmx_mtop_t *top_global, t_fcdata *fcd,
2360 t_state *state_global,
2362 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2363 gmx_edsam_t gmx_unused ed,
2365 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2366 gmx_membed_t gmx_unused membed,
2367 real gmx_unused cpt_period, real gmx_unused max_hours,
2368 const char gmx_unused *deviceOptions,
2369 unsigned long gmx_unused Flags,
2370 gmx_walltime_accounting_t walltime_accounting)
2372 const char *SD = "Steepest Descents";
2373 em_state_t *s_min, *s_try;
2375 gmx_localtop_t *top;
2376 gmx_enerdata_t *enerd;
2378 gmx_global_stat_t gstat;
2380 real stepsize, constepsize;
2384 gmx_bool bDone, bAbort, do_x, do_f;
2389 int steps_accepted = 0;
2393 s_min = init_em_state();
2394 s_try = init_em_state();
2396 /* Init em and store the local state in s_try */
2397 init_em(fplog, SD, cr, inputrec,
2398 state_global, top_global, s_try, &top, &f, &f_global,
2399 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2400 nfile, fnm, &outf, &mdebin);
2402 /* Print to log file */
2403 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2405 /* Set variables for stepsize (in nm). This is the largest
2406 * step that we are going to make in any direction.
2408 ustep = inputrec->em_stepsize;
2411 /* Max number of steps */
2412 nsteps = inputrec->nsteps;
2416 /* Print to the screen */
2417 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2421 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2424 /**** HERE STARTS THE LOOP ****
2425 * count is the counter for the number of steps
2426 * bDone will be TRUE when the minimization has converged
2427 * bAbort will be TRUE when nsteps steps have been performed or when
2428 * the stepsize becomes smaller than is reasonable for machine precision
2433 while (!bDone && !bAbort)
2435 bAbort = (nsteps >= 0) && (count == nsteps);
2437 /* set new coordinates, except for first step */
2440 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2441 s_min, stepsize, s_min->f, s_try,
2442 constr, top, nrnb, wcycle, count);
2445 evaluate_energy(fplog, cr,
2446 top_global, s_try, top,
2447 inputrec, nrnb, wcycle, gstat,
2448 vsite, constr, fcd, graph, mdatoms, fr,
2449 mu_tot, enerd, vir, pres, count, count == 0);
2453 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2458 s_min->epot = s_try->epot + 1;
2461 /* Print it if necessary */
2466 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2467 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2468 (s_try->epot < s_min->epot) ? '\n' : '\r');
2471 if (s_try->epot < s_min->epot)
2473 /* Store the new (lower) energies */
2474 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2475 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2476 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2477 print_ebin(outf->fp_ene, TRUE,
2478 do_per_step(steps_accepted, inputrec->nstdisreout),
2479 do_per_step(steps_accepted, inputrec->nstorireout),
2480 fplog, count, count, eprNORMAL, TRUE,
2481 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2486 /* Now if the new energy is smaller than the previous...
2487 * or if this is the first step!
2488 * or if we did random steps!
2491 if ( (count == 0) || (s_try->epot < s_min->epot) )
2495 /* Test whether the convergence criterion is met... */
2496 bDone = (s_try->fmax < inputrec->em_tol);
2498 /* Copy the arrays for force, positions and energy */
2499 /* The 'Min' array always holds the coords and forces of the minimal
2501 swap_em_state(s_min, s_try);
2507 /* Write to trn, if necessary */
2508 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2509 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2510 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2511 top_global, inputrec, count,
2512 s_min, state_global, f_global);
2516 /* If energy is not smaller make the step smaller... */
2519 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2521 /* Reload the old state */
2522 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2523 s_min, top, mdatoms, fr, vsite, constr,
2528 /* Determine new step */
2529 stepsize = ustep/s_min->fmax;
2531 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2533 if (count == nsteps || ustep < 1e-12)
2535 if (count == nsteps || ustep < 1e-6)
2540 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2541 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2547 } /* End of the loop */
2549 /* Print some shit... */
2552 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2554 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2555 top_global, inputrec, count,
2556 s_min, state_global, f_global);
2558 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2562 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2563 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2564 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2565 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2568 finish_em(cr, outf, walltime_accounting, wcycle);
2570 /* To print the actual number of steps we needed somewhere */
2571 inputrec->nsteps = count;
2573 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2576 } /* That's all folks */
2579 double do_nm(FILE *fplog, t_commrec *cr,
2580 int nfile, const t_filenm fnm[],
2581 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2582 int gmx_unused nstglobalcomm,
2583 gmx_vsite_t *vsite, gmx_constr_t constr,
2584 int gmx_unused stepout,
2585 t_inputrec *inputrec,
2586 gmx_mtop_t *top_global, t_fcdata *fcd,
2587 t_state *state_global,
2589 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2590 gmx_edsam_t gmx_unused ed,
2592 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2593 gmx_membed_t gmx_unused membed,
2594 real gmx_unused cpt_period, real gmx_unused max_hours,
2595 const char gmx_unused *deviceOptions,
2596 unsigned long gmx_unused Flags,
2597 gmx_walltime_accounting_t walltime_accounting)
2599 const char *NM = "Normal Mode Analysis";
2601 int natoms, atom, d;
2604 gmx_localtop_t *top;
2605 gmx_enerdata_t *enerd;
2607 gmx_global_stat_t gstat;
2609 real t, t0, lambda, lam0;
2614 gmx_bool bSparse; /* use sparse matrix storage format */
2616 gmx_sparsematrix_t * sparse_matrix = NULL;
2617 real * full_matrix = NULL;
2618 em_state_t * state_work;
2620 /* added with respect to mdrun */
2621 int i, j, k, row, col;
2622 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2628 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2631 state_work = init_em_state();
2633 /* Init em and store the local state in state_minimum */
2634 init_em(fplog, NM, cr, inputrec,
2635 state_global, top_global, state_work, &top,
2637 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2638 nfile, fnm, &outf, NULL);
2640 natoms = top_global->natoms;
2648 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2649 " which MIGHT not be accurate enough for normal mode analysis.\n"
2650 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2651 " are fairly modest even if you recompile in double precision.\n\n");
2655 /* Check if we can/should use sparse storage format.
2657 * Sparse format is only useful when the Hessian itself is sparse, which it
2658 * will be when we use a cutoff.
2659 * For small systems (n<1000) it is easier to always use full matrix format, though.
2661 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2663 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2666 else if (top_global->natoms < 1000)
2668 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2673 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2679 sz = DIM*top_global->natoms;
2681 fprintf(stderr, "Allocating Hessian memory...\n\n");
2685 sparse_matrix = gmx_sparsematrix_init(sz);
2686 sparse_matrix->compressed_symmetric = TRUE;
2690 snew(full_matrix, sz*sz);
2694 /* Initial values */
2695 t0 = inputrec->init_t;
2696 lam0 = inputrec->fepvals->init_lambda;
2704 /* Write start time and temperature */
2705 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2707 /* fudge nr of steps to nr of atoms */
2708 inputrec->nsteps = natoms*2;
2712 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2713 *(top_global->name), (int)inputrec->nsteps);
2716 nnodes = cr->nnodes;
2718 /* Make evaluate_energy do a single node force calculation */
2720 evaluate_energy(fplog, cr,
2721 top_global, state_work, top,
2722 inputrec, nrnb, wcycle, gstat,
2723 vsite, constr, fcd, graph, mdatoms, fr,
2724 mu_tot, enerd, vir, pres, -1, TRUE);
2725 cr->nnodes = nnodes;
2727 /* if forces are not small, warn user */
2728 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2730 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2731 if (state_work->fmax > 1.0e-3)
2733 md_print_info(cr, fplog,
2734 "The force is probably not small enough to "
2735 "ensure that you are at a minimum.\n"
2736 "Be aware that negative eigenvalues may occur\n"
2737 "when the resulting matrix is diagonalized.\n\n");
2740 /***********************************************************
2742 * Loop over all pairs in matrix
2744 * do_force called twice. Once with positive and
2745 * once with negative displacement
2747 ************************************************************/
2749 /* Steps are divided one by one over the nodes */
2750 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2753 for (d = 0; d < DIM; d++)
2755 x_min = state_work->s.x[atom][d];
2757 state_work->s.x[atom][d] = x_min - der_range;
2759 /* Make evaluate_energy do a single node force calculation */
2761 evaluate_energy(fplog, cr,
2762 top_global, state_work, top,
2763 inputrec, nrnb, wcycle, gstat,
2764 vsite, constr, fcd, graph, mdatoms, fr,
2765 mu_tot, enerd, vir, pres, atom*2, FALSE);
2767 for (i = 0; i < natoms; i++)
2769 copy_rvec(state_work->f[i], fneg[i]);
2772 state_work->s.x[atom][d] = x_min + der_range;
2774 evaluate_energy(fplog, cr,
2775 top_global, state_work, top,
2776 inputrec, nrnb, wcycle, gstat,
2777 vsite, constr, fcd, graph, mdatoms, fr,
2778 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2779 cr->nnodes = nnodes;
2781 /* x is restored to original */
2782 state_work->s.x[atom][d] = x_min;
2784 for (j = 0; j < natoms; j++)
2786 for (k = 0; (k < DIM); k++)
2789 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2797 #define mpi_type MPI_DOUBLE
2799 #define mpi_type MPI_FLOAT
2801 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2802 cr->mpi_comm_mygroup);
2807 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2813 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2814 cr->mpi_comm_mygroup, &stat);
2819 row = (atom + node)*DIM + d;
2821 for (j = 0; j < natoms; j++)
2823 for (k = 0; k < DIM; k++)
2829 if (col >= row && dfdx[j][k] != 0.0)
2831 gmx_sparsematrix_increment_value(sparse_matrix,
2832 row, col, dfdx[j][k]);
2837 full_matrix[row*sz+col] = dfdx[j][k];
2844 if (bVerbose && fplog)
2849 /* write progress */
2850 if (MASTER(cr) && bVerbose)
2852 fprintf(stderr, "\rFinished step %d out of %d",
2853 min(atom+nnodes, natoms), natoms);
2860 fprintf(stderr, "\n\nWriting Hessian...\n");
2861 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2864 finish_em(cr, outf, walltime_accounting, wcycle);
2866 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);