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
46 #include "gromacs/fileio/confio.h"
54 #include "gmx_fatal.h"
66 #include "md_support.h"
69 #include "gromacs/fileio/trnio.h"
72 #include "gmx_wallcycle.h"
73 #include "mtop_util.h"
74 #include "gromacs/fileio/gmxfio.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/fileio/trajectory_writing.h"
84 #include "gromacs/timing/walltime_accounting.h"
95 static em_state_t *init_em_state()
101 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
102 snew(ems->s.lambda, efptNR);
107 static void print_em_start(FILE *fplog,
109 gmx_walltime_accounting_t walltime_accounting,
110 gmx_wallcycle_t wcycle,
115 walltime_accounting_start(walltime_accounting);
117 sprintf(buf, "Started %s", name);
118 print_date_and_time(fplog, cr->nodeid, buf, NULL);
120 wallcycle_start(wcycle, ewcRUN);
122 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
123 gmx_wallcycle_t wcycle)
125 wallcycle_stop(wcycle, ewcRUN);
127 walltime_accounting_end(walltime_accounting);
130 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
133 fprintf(out, "%s:\n", minimizer);
134 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
135 fprintf(out, " Number of steps = %12d\n", nsteps);
138 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
144 "\nEnergy minimization reached the maximum number "
145 "of steps before the forces reached the requested "
146 "precision Fmax < %g.\n", ftol);
151 "\nEnergy minimization has stopped, but the forces have "
152 "not converged to the requested precision Fmax < %g (which "
153 "may not be possible for your system). It stopped "
154 "because the algorithm tried to make a new step whose size "
155 "was too small, or there was no change in the energy since "
156 "last step. Either way, we regard the minimization as "
157 "converged to within the available machine precision, "
158 "given your starting configuration and EM parameters.\n%s%s",
160 sizeof(real) < sizeof(double) ?
161 "\nDouble precision normally gives you higher accuracy, but "
162 "this is often not needed for preparing to run molecular "
166 "You might need to increase your constraint accuracy, or turn\n"
167 "off constraints altogether (set constraints = none in mdp file)\n" :
170 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
175 static void print_converged(FILE *fp, const char *alg, real ftol,
176 gmx_large_int_t count, gmx_bool bDone, gmx_large_int_t nsteps,
177 real epot, real fmax, int nfmax, real fnorm)
179 char buf[STEPSTRSIZE];
183 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
184 alg, ftol, gmx_step_str(count, buf));
186 else if (count < nsteps)
188 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
189 "but did not reach the requested Fmax < %g.\n",
190 alg, gmx_step_str(count, buf), ftol);
194 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
195 alg, ftol, gmx_step_str(count, buf));
199 fprintf(fp, "Potential Energy = %21.14e\n", epot);
200 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
201 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
203 fprintf(fp, "Potential Energy = %14.7e\n", epot);
204 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
205 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
209 static void get_f_norm_max(t_commrec *cr,
210 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
211 real *fnorm, real *fmax, int *a_fmax)
214 real fmax2, fmax2_0, fam;
215 int la_max, a_max, start, end, i, m, gf;
217 /* This routine finds the largest force and returns it.
218 * On parallel machines the global max is taken.
224 start = mdatoms->start;
225 end = mdatoms->homenr + start;
226 if (mdatoms->cFREEZE)
228 for (i = start; i < end; i++)
230 gf = mdatoms->cFREEZE[i];
232 for (m = 0; m < DIM; m++)
234 if (!opts->nFreeze[gf][m])
249 for (i = start; i < end; i++)
261 if (la_max >= 0 && DOMAINDECOMP(cr))
263 a_max = cr->dd->gatindex[la_max];
271 snew(sum, 2*cr->nnodes+1);
272 sum[2*cr->nodeid] = fmax2;
273 sum[2*cr->nodeid+1] = a_max;
274 sum[2*cr->nnodes] = fnorm2;
275 gmx_sumd(2*cr->nnodes+1, sum, cr);
276 fnorm2 = sum[2*cr->nnodes];
277 /* Determine the global maximum */
278 for (i = 0; i < cr->nnodes; i++)
280 if (sum[2*i] > fmax2)
283 a_max = (int)(sum[2*i+1] + 0.5);
291 *fnorm = sqrt(fnorm2);
303 static void get_state_f_norm_max(t_commrec *cr,
304 t_grpopts *opts, t_mdatoms *mdatoms,
307 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
310 void init_em(FILE *fplog, const char *title,
311 t_commrec *cr, t_inputrec *ir,
312 t_state *state_global, gmx_mtop_t *top_global,
313 em_state_t *ems, gmx_localtop_t **top,
314 rvec **f, rvec **f_global,
315 t_nrnb *nrnb, rvec mu_tot,
316 t_forcerec *fr, gmx_enerdata_t **enerd,
317 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
318 gmx_vsite_t *vsite, gmx_constr_t constr,
319 int nfile, const t_filenm fnm[],
320 gmx_mdoutf_t **outf, t_mdebin **mdebin)
322 int start, homenr, i;
327 fprintf(fplog, "Initiating %s\n", title);
330 state_global->ngtc = 0;
332 /* Initialize lambda variables */
333 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
337 if (DOMAINDECOMP(cr))
339 *top = dd_init_local_top(top_global);
341 dd_init_local_state(cr->dd, state_global, &ems->s);
345 /* Distribute the charge groups over the nodes from the master node */
346 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
347 state_global, top_global, ir,
348 &ems->s, &ems->f, mdatoms, *top,
349 fr, vsite, NULL, constr,
351 dd_store_state(cr->dd, &ems->s);
355 snew(*f_global, top_global->natoms);
365 snew(*f, top_global->natoms);
367 /* Just copy the state */
368 ems->s = *state_global;
369 snew(ems->s.x, ems->s.nalloc);
370 snew(ems->f, ems->s.nalloc);
371 for (i = 0; i < state_global->natoms; i++)
373 copy_rvec(state_global->x[i], ems->s.x[i]);
375 copy_mat(state_global->box, ems->s.box);
377 if (PAR(cr) && ir->eI != eiNM)
379 /* Initialize the particle decomposition and split the topology */
380 *top = split_system(fplog, top_global, ir, cr);
382 pd_cg_range(cr, &fr->cg0, &fr->hcg);
386 *top = gmx_mtop_generate_local_top(top_global, ir);
390 forcerec_set_excl_load(fr, *top, cr);
392 setup_bonded_threading(fr, &(*top)->idef);
394 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
396 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
405 pd_at_range(cr, &start, &homenr);
411 homenr = top_global->natoms;
413 atoms2md(top_global, ir, 0, NULL, start, homenr, mdatoms);
414 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
418 set_vsite_top(vsite, *top, mdatoms, cr);
424 if (ir->eConstrAlg == econtSHAKE &&
425 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
427 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
428 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
431 if (!DOMAINDECOMP(cr))
433 set_constraints(constr, *top, ir, mdatoms, cr);
436 if (!ir->bContinuation)
438 /* Constrain the starting coordinates */
440 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
441 ir, NULL, cr, -1, 0, mdatoms,
442 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
443 ems->s.lambda[efptFEP], &dvdl_constr,
444 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
450 *gstat = global_stat_init(ir);
453 *outf = init_mdoutf(nfile, fnm, 0, cr, ir, NULL);
456 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
461 /* Init bin for energy stuff */
462 *mdebin = init_mdebin((*outf)->fp_ene, top_global, ir, NULL);
466 calc_shifts(ems->s.box, fr->shift_vec);
469 static void finish_em(t_commrec *cr, gmx_mdoutf_t *outf,
470 gmx_walltime_accounting_t walltime_accounting,
471 gmx_wallcycle_t wcycle)
473 if (!(cr->duty & DUTY_PME))
475 /* Tell the PME only node to finish */
476 gmx_pme_send_finish(cr);
481 em_time_end(walltime_accounting, wcycle);
484 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
493 static void copy_em_coords(em_state_t *ems, t_state *state)
497 for (i = 0; (i < state->natoms); i++)
499 copy_rvec(ems->s.x[i], state->x[i]);
503 static void write_em_traj(FILE *fplog, t_commrec *cr,
505 gmx_bool bX, gmx_bool bF, const char *confout,
506 gmx_mtop_t *top_global,
507 t_inputrec *ir, gmx_large_int_t step,
509 t_state *state_global, rvec *f_global)
513 if ((bX || bF || confout != NULL) && !DOMAINDECOMP(cr))
515 copy_em_coords(state, state_global);
522 mdof_flags |= MDOF_X;
526 mdof_flags |= MDOF_F;
528 write_traj(fplog, cr, outf, mdof_flags,
529 top_global, step, (double)step,
530 &state->s, state_global, state->f, f_global, NULL, NULL);
532 if (confout != NULL && MASTER(cr))
534 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
536 /* Make molecules whole only for confout writing */
537 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
541 write_sto_conf_mtop(confout,
542 *top_global->name, top_global,
543 state_global->x, NULL, ir->ePBC, state_global->box);
547 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
549 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
550 gmx_constr_t constr, gmx_localtop_t *top,
551 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
552 gmx_large_int_t count)
564 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
566 gmx_incons("state mismatch in do_em_step");
569 s2->flags = s1->flags;
571 if (s2->nalloc != s1->nalloc)
573 s2->nalloc = s1->nalloc;
574 srenew(s2->x, s1->nalloc);
575 srenew(ems2->f, s1->nalloc);
576 if (s2->flags & (1<<estCGP))
578 srenew(s2->cg_p, s1->nalloc);
582 s2->natoms = s1->natoms;
583 copy_mat(s1->box, s2->box);
584 /* Copy free energy state */
585 for (i = 0; i < efptNR; i++)
587 s2->lambda[i] = s1->lambda[i];
589 copy_mat(s1->box, s2->box);
592 end = md->start + md->homenr;
597 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
602 #pragma omp for schedule(static) nowait
603 for (i = start; i < end; i++)
609 for (m = 0; m < DIM; m++)
611 if (ir->opts.nFreeze[gf][m])
617 x2[i][m] = x1[i][m] + a*f[i][m];
622 if (s2->flags & (1<<estCGP))
624 /* Copy the CG p vector */
627 #pragma omp for schedule(static) nowait
628 for (i = start; i < end; i++)
630 copy_rvec(x1[i], x2[i]);
634 if (DOMAINDECOMP(cr))
636 s2->ddp_count = s1->ddp_count;
637 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
640 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
641 srenew(s2->cg_gl, s2->cg_gl_nalloc);
644 s2->ncg_gl = s1->ncg_gl;
645 #pragma omp for schedule(static) nowait
646 for (i = 0; i < s2->ncg_gl; i++)
648 s2->cg_gl[i] = s1->cg_gl[i];
650 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
656 wallcycle_start(wcycle, ewcCONSTR);
658 constrain(NULL, TRUE, TRUE, constr, &top->idef,
659 ir, NULL, cr, count, 0, md,
660 s1->x, s2->x, NULL, bMolPBC, s2->box,
661 s2->lambda[efptBONDED], &dvdl_constr,
662 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
663 wallcycle_stop(wcycle, ewcCONSTR);
667 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
668 gmx_mtop_t *top_global, t_inputrec *ir,
669 em_state_t *ems, gmx_localtop_t *top,
670 t_mdatoms *mdatoms, t_forcerec *fr,
671 gmx_vsite_t *vsite, gmx_constr_t constr,
672 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
674 /* Repartition the domain decomposition */
675 wallcycle_start(wcycle, ewcDOMDEC);
676 dd_partition_system(fplog, step, cr, FALSE, 1,
677 NULL, top_global, ir,
679 mdatoms, top, fr, vsite, NULL, constr,
680 nrnb, wcycle, FALSE);
681 dd_store_state(cr->dd, &ems->s);
682 wallcycle_stop(wcycle, ewcDOMDEC);
685 static void evaluate_energy(FILE *fplog, t_commrec *cr,
686 gmx_mtop_t *top_global,
687 em_state_t *ems, gmx_localtop_t *top,
688 t_inputrec *inputrec,
689 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
690 gmx_global_stat_t gstat,
691 gmx_vsite_t *vsite, gmx_constr_t constr,
693 t_graph *graph, t_mdatoms *mdatoms,
694 t_forcerec *fr, rvec mu_tot,
695 gmx_enerdata_t *enerd, tensor vir, tensor pres,
696 gmx_large_int_t count, gmx_bool bFirst)
701 tensor force_vir, shake_vir, ekin;
702 real dvdl_constr, prescorr, enercorr, dvdlcorr;
705 /* Set the time to the initial time, the time does not change during EM */
706 t = inputrec->init_t;
709 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
711 /* This the first state or an old state used before the last ns */
717 if (inputrec->nstlist > 0)
721 else if (inputrec->nstlist == -1)
723 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
726 gmx_sumi(1, &nabnsb, cr);
734 construct_vsites(vsite, ems->s.x, 1, NULL,
735 top->idef.iparams, top->idef.il,
736 fr->ePBC, fr->bMolPBC, graph, cr, ems->s.box);
739 if (DOMAINDECOMP(cr))
743 /* Repartition the domain decomposition */
744 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
745 ems, top, mdatoms, fr, vsite, constr,
750 /* Calc force & energy on new trial position */
751 /* do_force always puts the charge groups in the box and shifts again
752 * We do not unshift, so molecules are always whole in congrad.c
754 do_force(fplog, cr, inputrec,
755 count, nrnb, wcycle, top, &top_global->groups,
756 ems->s.box, ems->s.x, &ems->s.hist,
757 ems->f, force_vir, mdatoms, enerd, fcd,
758 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
759 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
760 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
761 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
763 /* Clear the unused shake virial and pressure */
764 clear_mat(shake_vir);
767 /* Communicate stuff when parallel */
768 if (PAR(cr) && inputrec->eI != eiNM)
770 wallcycle_start(wcycle, ewcMoveE);
772 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
773 inputrec, NULL, NULL, NULL, 1, &terminate,
774 top_global, &ems->s, FALSE,
780 wallcycle_stop(wcycle, ewcMoveE);
783 /* Calculate long range corrections to pressure and energy */
784 calc_dispcorr(fplog, inputrec, fr, count, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
785 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
786 enerd->term[F_DISPCORR] = enercorr;
787 enerd->term[F_EPOT] += enercorr;
788 enerd->term[F_PRES] += prescorr;
789 enerd->term[F_DVDL] += dvdlcorr;
791 ems->epot = enerd->term[F_EPOT];
795 /* Project out the constraint components of the force */
796 wallcycle_start(wcycle, ewcCONSTR);
798 constrain(NULL, FALSE, FALSE, constr, &top->idef,
799 inputrec, NULL, cr, count, 0, mdatoms,
800 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
801 ems->s.lambda[efptBONDED], &dvdl_constr,
802 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
803 if (fr->bSepDVDL && fplog)
805 gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
807 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
808 m_add(force_vir, shake_vir, vir);
809 wallcycle_stop(wcycle, ewcCONSTR);
813 copy_mat(force_vir, vir);
817 enerd->term[F_PRES] =
818 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
820 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
822 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
824 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
828 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
830 em_state_t *s_min, em_state_t *s_b)
834 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
836 unsigned char *grpnrFREEZE;
840 fprintf(debug, "Doing reorder_partsum\n");
846 cgs_gl = dd_charge_groups_global(cr->dd);
847 index = cgs_gl->index;
849 /* Collect fm in a global vector fmg.
850 * This conflicts with the spirit of domain decomposition,
851 * but to fully optimize this a much more complicated algorithm is required.
853 snew(fmg, mtop->natoms);
855 ncg = s_min->s.ncg_gl;
856 cg_gl = s_min->s.cg_gl;
858 for (c = 0; c < ncg; c++)
863 for (a = a0; a < a1; a++)
865 copy_rvec(fm[i], fmg[a]);
869 gmx_sum(mtop->natoms*3, fmg[0], cr);
871 /* Now we will determine the part of the sum for the cgs in state s_b */
873 cg_gl = s_b->s.cg_gl;
877 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
878 for (c = 0; c < ncg; c++)
883 for (a = a0; a < a1; a++)
885 if (mdatoms->cFREEZE && grpnrFREEZE)
889 for (m = 0; m < DIM; m++)
891 if (!opts->nFreeze[gf][m])
893 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
905 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
907 em_state_t *s_min, em_state_t *s_b)
913 /* This is just the classical Polak-Ribiere calculation of beta;
914 * it looks a bit complicated since we take freeze groups into account,
915 * and might have to sum it in parallel runs.
918 if (!DOMAINDECOMP(cr) ||
919 (s_min->s.ddp_count == cr->dd->ddp_count &&
920 s_b->s.ddp_count == cr->dd->ddp_count))
926 /* This part of code can be incorrect with DD,
927 * since the atom ordering in s_b and s_min might differ.
929 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
931 if (mdatoms->cFREEZE)
933 gf = mdatoms->cFREEZE[i];
935 for (m = 0; m < DIM; m++)
937 if (!opts->nFreeze[gf][m])
939 sum += (fb[i][m] - fm[i][m])*fb[i][m];
946 /* We need to reorder cgs while summing */
947 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
951 gmx_sumd(1, &sum, cr);
954 return sum/sqr(s_min->fnorm);
957 double do_cg(FILE *fplog, t_commrec *cr,
958 int nfile, const t_filenm fnm[],
959 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
960 int gmx_unused nstglobalcomm,
961 gmx_vsite_t *vsite, gmx_constr_t constr,
962 int gmx_unused stepout,
963 t_inputrec *inputrec,
964 gmx_mtop_t *top_global, t_fcdata *fcd,
965 t_state *state_global,
967 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
968 gmx_edsam_t gmx_unused ed,
970 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
971 gmx_membed_t gmx_unused membed,
972 real gmx_unused cpt_period, real gmx_unused max_hours,
973 const char gmx_unused *deviceOptions,
974 unsigned long gmx_unused Flags,
975 gmx_walltime_accounting_t walltime_accounting)
977 const char *CG = "Polak-Ribiere Conjugate Gradients";
979 em_state_t *s_min, *s_a, *s_b, *s_c;
981 gmx_enerdata_t *enerd;
983 gmx_global_stat_t gstat;
985 rvec *f_global, *p, *sf, *sfm;
986 double gpa, gpb, gpc, tmp, sum[2], minstep;
989 real a, b, c, beta = 0.0;
993 gmx_bool converged, foundlower;
995 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
997 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
999 int i, m, gf, step, nminstep;
1004 s_min = init_em_state();
1005 s_a = init_em_state();
1006 s_b = init_em_state();
1007 s_c = init_em_state();
1009 /* Init em and store the local state in s_min */
1010 init_em(fplog, CG, cr, inputrec,
1011 state_global, top_global, s_min, &top, &f, &f_global,
1012 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1013 nfile, fnm, &outf, &mdebin);
1015 /* Print to log file */
1016 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1018 /* Max number of steps */
1019 number_steps = inputrec->nsteps;
1023 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1027 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1030 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1031 /* do_force always puts the charge groups in the box and shifts again
1032 * We do not unshift, so molecules are always whole in congrad.c
1034 evaluate_energy(fplog, cr,
1035 top_global, s_min, top,
1036 inputrec, nrnb, wcycle, gstat,
1037 vsite, constr, fcd, graph, mdatoms, fr,
1038 mu_tot, enerd, vir, pres, -1, TRUE);
1043 /* Copy stuff to the energy bin for easy printing etc. */
1044 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1045 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1046 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1048 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1049 print_ebin(outf->fp_ene, TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1050 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1054 /* Estimate/guess the initial stepsize */
1055 stepsize = inputrec->em_stepsize/s_min->fnorm;
1059 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1060 s_min->fmax, s_min->a_fmax+1);
1061 fprintf(stderr, " F-Norm = %12.5e\n",
1062 s_min->fnorm/sqrt(state_global->natoms));
1063 fprintf(stderr, "\n");
1064 /* and copy to the log file too... */
1065 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1066 s_min->fmax, s_min->a_fmax+1);
1067 fprintf(fplog, " F-Norm = %12.5e\n",
1068 s_min->fnorm/sqrt(state_global->natoms));
1069 fprintf(fplog, "\n");
1071 /* Start the loop over CG steps.
1072 * Each successful step is counted, and we continue until
1073 * we either converge or reach the max number of steps.
1076 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1079 /* start taking steps in a new direction
1080 * First time we enter the routine, beta=0, and the direction is
1081 * simply the negative gradient.
1084 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1089 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1091 if (mdatoms->cFREEZE)
1093 gf = mdatoms->cFREEZE[i];
1095 for (m = 0; m < DIM; m++)
1097 if (!inputrec->opts.nFreeze[gf][m])
1099 p[i][m] = sf[i][m] + beta*p[i][m];
1100 gpa -= p[i][m]*sf[i][m];
1101 /* f is negative gradient, thus the sign */
1110 /* Sum the gradient along the line across CPUs */
1113 gmx_sumd(1, &gpa, cr);
1116 /* Calculate the norm of the search vector */
1117 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1119 /* Just in case stepsize reaches zero due to numerical precision... */
1122 stepsize = inputrec->em_stepsize/pnorm;
1126 * Double check the value of the derivative in the search direction.
1127 * If it is positive it must be due to the old information in the
1128 * CG formula, so just remove that and start over with beta=0.
1129 * This corresponds to a steepest descent step.
1134 step--; /* Don't count this step since we are restarting */
1135 continue; /* Go back to the beginning of the big for-loop */
1138 /* Calculate minimum allowed stepsize, before the average (norm)
1139 * relative change in coordinate is smaller than precision
1142 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1144 for (m = 0; m < DIM; m++)
1146 tmp = fabs(s_min->s.x[i][m]);
1155 /* Add up from all CPUs */
1158 gmx_sumd(1, &minstep, cr);
1161 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1163 if (stepsize < minstep)
1169 /* Write coordinates if necessary */
1170 do_x = do_per_step(step, inputrec->nstxout);
1171 do_f = do_per_step(step, inputrec->nstfout);
1173 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1174 top_global, inputrec, step,
1175 s_min, state_global, f_global);
1177 /* Take a step downhill.
1178 * In theory, we should minimize the function along this direction.
1179 * That is quite possible, but it turns out to take 5-10 function evaluations
1180 * for each line. However, we dont really need to find the exact minimum -
1181 * it is much better to start a new CG step in a modified direction as soon
1182 * as we are close to it. This will save a lot of energy evaluations.
1184 * In practice, we just try to take a single step.
1185 * If it worked (i.e. lowered the energy), we increase the stepsize but
1186 * the continue straight to the next CG step without trying to find any minimum.
1187 * If it didn't work (higher energy), there must be a minimum somewhere between
1188 * the old position and the new one.
1190 * Due to the finite numerical accuracy, it turns out that it is a good idea
1191 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1192 * This leads to lower final energies in the tests I've done. / Erik
1194 s_a->epot = s_min->epot;
1196 c = a + stepsize; /* reference position along line is zero */
1198 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1200 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1201 s_min, top, mdatoms, fr, vsite, constr,
1205 /* Take a trial step (new coords in s_c) */
1206 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1207 constr, top, nrnb, wcycle, -1);
1210 /* Calculate energy for the trial step */
1211 evaluate_energy(fplog, cr,
1212 top_global, s_c, top,
1213 inputrec, nrnb, wcycle, gstat,
1214 vsite, constr, fcd, graph, mdatoms, fr,
1215 mu_tot, enerd, vir, pres, -1, FALSE);
1217 /* Calc derivative along line */
1221 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1223 for (m = 0; m < DIM; m++)
1225 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1228 /* Sum the gradient along the line across CPUs */
1231 gmx_sumd(1, &gpc, cr);
1234 /* This is the max amount of increase in energy we tolerate */
1235 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1237 /* Accept the step if the energy is lower, or if it is not significantly higher
1238 * and the line derivative is still negative.
1240 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1243 /* Great, we found a better energy. Increase step for next iteration
1244 * if we are still going down, decrease it otherwise
1248 stepsize *= 1.618034; /* The golden section */
1252 stepsize *= 0.618034; /* 1/golden section */
1257 /* New energy is the same or higher. We will have to do some work
1258 * to find a smaller value in the interval. Take smaller step next time!
1261 stepsize *= 0.618034;
1267 /* OK, if we didn't find a lower value we will have to locate one now - there must
1268 * be one in the interval [a=0,c].
1269 * The same thing is valid here, though: Don't spend dozens of iterations to find
1270 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1271 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1273 * I also have a safeguard for potentially really patological functions so we never
1274 * take more than 20 steps before we give up ...
1276 * If we already found a lower value we just skip this step and continue to the update.
1284 /* Select a new trial point.
1285 * If the derivatives at points a & c have different sign we interpolate to zero,
1286 * otherwise just do a bisection.
1288 if (gpa < 0 && gpc > 0)
1290 b = a + gpa*(a-c)/(gpc-gpa);
1297 /* safeguard if interpolation close to machine accuracy causes errors:
1298 * never go outside the interval
1300 if (b <= a || b >= c)
1305 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1307 /* Reload the old state */
1308 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1309 s_min, top, mdatoms, fr, vsite, constr,
1313 /* Take a trial step to this new point - new coords in s_b */
1314 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1315 constr, top, nrnb, wcycle, -1);
1318 /* Calculate energy for the trial step */
1319 evaluate_energy(fplog, cr,
1320 top_global, s_b, top,
1321 inputrec, nrnb, wcycle, gstat,
1322 vsite, constr, fcd, graph, mdatoms, fr,
1323 mu_tot, enerd, vir, pres, -1, FALSE);
1325 /* p does not change within a step, but since the domain decomposition
1326 * might change, we have to use cg_p of s_b here.
1331 for (i = mdatoms->start; i < mdatoms->start+mdatoms->homenr; i++)
1333 for (m = 0; m < DIM; m++)
1335 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1338 /* Sum the gradient along the line across CPUs */
1341 gmx_sumd(1, &gpb, cr);
1346 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1347 s_a->epot, s_b->epot, s_c->epot, gpb);
1350 epot_repl = s_b->epot;
1352 /* Keep one of the intervals based on the value of the derivative at the new point */
1355 /* Replace c endpoint with b */
1356 swap_em_state(s_b, s_c);
1362 /* Replace a endpoint with b */
1363 swap_em_state(s_b, s_a);
1369 * Stop search as soon as we find a value smaller than the endpoints.
1370 * Never run more than 20 steps, no matter what.
1374 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1377 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1380 /* OK. We couldn't find a significantly lower energy.
1381 * If beta==0 this was steepest descent, and then we give up.
1382 * If not, set beta=0 and restart with steepest descent before quitting.
1392 /* Reset memory before giving up */
1398 /* Select min energy state of A & C, put the best in B.
1400 if (s_c->epot < s_a->epot)
1404 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1405 s_c->epot, s_a->epot);
1407 swap_em_state(s_b, s_c);
1415 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1416 s_a->epot, s_c->epot);
1418 swap_em_state(s_b, s_a);
1428 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1431 swap_em_state(s_b, s_c);
1436 /* new search direction */
1437 /* beta = 0 means forget all memory and restart with steepest descents. */
1438 if (nstcg && ((step % nstcg) == 0))
1444 /* s_min->fnorm cannot be zero, because then we would have converged
1448 /* Polak-Ribiere update.
1449 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1451 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1453 /* Limit beta to prevent oscillations */
1454 if (fabs(beta) > 5.0)
1460 /* update positions */
1461 swap_em_state(s_min, s_b);
1464 /* Print it if necessary */
1469 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1470 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1471 s_min->fmax, s_min->a_fmax+1);
1473 /* Store the new (lower) energies */
1474 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1475 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1476 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1478 do_log = do_per_step(step, inputrec->nstlog);
1479 do_ene = do_per_step(step, inputrec->nstenergy);
1482 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1484 print_ebin(outf->fp_ene, do_ene, FALSE, FALSE,
1485 do_log ? fplog : NULL, step, step, eprNORMAL,
1486 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1489 /* Stop when the maximum force lies below tolerance.
1490 * If we have reached machine precision, converged is already set to true.
1492 converged = converged || (s_min->fmax < inputrec->em_tol);
1494 } /* End of the loop */
1498 step--; /* we never took that last step in this case */
1501 if (s_min->fmax > inputrec->em_tol)
1505 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1506 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1513 /* If we printed energy and/or logfile last step (which was the last step)
1514 * we don't have to do it again, but otherwise print the final values.
1518 /* Write final value to log since we didn't do anything the last step */
1519 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1521 if (!do_ene || !do_log)
1523 /* Write final energy file entries */
1524 print_ebin(outf->fp_ene, !do_ene, FALSE, FALSE,
1525 !do_log ? fplog : NULL, step, step, eprNORMAL,
1526 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1530 /* Print some stuff... */
1533 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1537 * For accurate normal mode calculation it is imperative that we
1538 * store the last conformation into the full precision binary trajectory.
1540 * However, we should only do it if we did NOT already write this step
1541 * above (which we did if do_x or do_f was true).
1543 do_x = !do_per_step(step, inputrec->nstxout);
1544 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1546 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1547 top_global, inputrec, step,
1548 s_min, state_global, f_global);
1550 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1554 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1555 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1556 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1557 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1559 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1562 finish_em(cr, outf, walltime_accounting, wcycle);
1564 /* To print the actual number of steps we needed somewhere */
1565 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1568 } /* That's all folks */
1571 double do_lbfgs(FILE *fplog, t_commrec *cr,
1572 int nfile, const t_filenm fnm[],
1573 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1574 int gmx_unused nstglobalcomm,
1575 gmx_vsite_t *vsite, gmx_constr_t constr,
1576 int gmx_unused stepout,
1577 t_inputrec *inputrec,
1578 gmx_mtop_t *top_global, t_fcdata *fcd,
1581 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1582 gmx_edsam_t gmx_unused ed,
1584 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1585 gmx_membed_t gmx_unused membed,
1586 real gmx_unused cpt_period, real gmx_unused max_hours,
1587 const char gmx_unused *deviceOptions,
1588 unsigned long gmx_unused Flags,
1589 gmx_walltime_accounting_t walltime_accounting)
1591 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1593 gmx_localtop_t *top;
1594 gmx_enerdata_t *enerd;
1596 gmx_global_stat_t gstat;
1599 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1600 double stepsize, gpa, gpb, gpc, tmp, minstep;
1601 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1602 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1603 real a, b, c, maxdelta, delta;
1604 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1605 real dgdx, dgdg, sq, yr, beta;
1607 gmx_bool converged, first;
1610 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1612 int start, end, number_steps;
1614 int i, k, m, n, nfmax, gf, step;
1621 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1626 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).");
1629 n = 3*state->natoms;
1630 nmaxcorr = inputrec->nbfgscorr;
1632 /* Allocate memory */
1633 /* Use pointers to real so we dont have to loop over both atoms and
1634 * dimensions all the time...
1635 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1636 * that point to the same memory.
1649 snew(rho, nmaxcorr);
1650 snew(alpha, nmaxcorr);
1653 for (i = 0; i < nmaxcorr; i++)
1659 for (i = 0; i < nmaxcorr; i++)
1668 init_em(fplog, LBFGS, cr, inputrec,
1669 state, top_global, &ems, &top, &f, &f_global,
1670 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1671 nfile, fnm, &outf, &mdebin);
1672 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1673 * so we free some memory again.
1678 xx = (real *)state->x;
1681 start = mdatoms->start;
1682 end = mdatoms->homenr + start;
1684 /* Print to log file */
1685 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1687 do_log = do_ene = do_x = do_f = TRUE;
1689 /* Max number of steps */
1690 number_steps = inputrec->nsteps;
1692 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1694 for (i = start; i < end; i++)
1696 if (mdatoms->cFREEZE)
1698 gf = mdatoms->cFREEZE[i];
1700 for (m = 0; m < DIM; m++)
1702 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1707 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1711 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1716 construct_vsites(vsite, state->x, 1, NULL,
1717 top->idef.iparams, top->idef.il,
1718 fr->ePBC, fr->bMolPBC, graph, cr, state->box);
1721 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1722 /* do_force always puts the charge groups in the box and shifts again
1723 * We do not unshift, so molecules are always whole
1728 evaluate_energy(fplog, cr,
1729 top_global, &ems, top,
1730 inputrec, nrnb, wcycle, gstat,
1731 vsite, constr, fcd, graph, mdatoms, fr,
1732 mu_tot, enerd, vir, pres, -1, TRUE);
1737 /* Copy stuff to the energy bin for easy printing etc. */
1738 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1739 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1740 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1742 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1743 print_ebin(outf->fp_ene, TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1744 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1748 /* This is the starting energy */
1749 Epot = enerd->term[F_EPOT];
1755 /* Set the initial step.
1756 * since it will be multiplied by the non-normalized search direction
1757 * vector (force vector the first time), we scale it by the
1758 * norm of the force.
1763 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1764 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1765 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1766 fprintf(stderr, "\n");
1767 /* and copy to the log file too... */
1768 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1769 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1770 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1771 fprintf(fplog, "\n");
1775 for (i = 0; i < n; i++)
1779 dx[point][i] = ff[i]; /* Initial search direction */
1787 stepsize = 1.0/fnorm;
1790 /* Start the loop over BFGS steps.
1791 * Each successful step is counted, and we continue until
1792 * we either converge or reach the max number of steps.
1797 /* Set the gradient from the force */
1799 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1802 /* Write coordinates if necessary */
1803 do_x = do_per_step(step, inputrec->nstxout);
1804 do_f = do_per_step(step, inputrec->nstfout);
1809 mdof_flags |= MDOF_X;
1814 mdof_flags |= MDOF_F;
1817 write_traj(fplog, cr, outf, mdof_flags,
1818 top_global, step, (real)step, state, state, f, f, NULL, NULL);
1820 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1822 /* pointer to current direction - point=0 first time here */
1825 /* calculate line gradient */
1826 for (gpa = 0, i = 0; i < n; i++)
1831 /* Calculate minimum allowed stepsize, before the average (norm)
1832 * relative change in coordinate is smaller than precision
1834 for (minstep = 0, i = 0; i < n; i++)
1844 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1846 if (stepsize < minstep)
1852 /* Store old forces and coordinates */
1853 for (i = 0; i < n; i++)
1862 for (i = 0; i < n; i++)
1867 /* Take a step downhill.
1868 * In theory, we should minimize the function along this direction.
1869 * That is quite possible, but it turns out to take 5-10 function evaluations
1870 * for each line. However, we dont really need to find the exact minimum -
1871 * it is much better to start a new BFGS step in a modified direction as soon
1872 * as we are close to it. This will save a lot of energy evaluations.
1874 * In practice, we just try to take a single step.
1875 * If it worked (i.e. lowered the energy), we increase the stepsize but
1876 * the continue straight to the next BFGS step without trying to find any minimum.
1877 * If it didn't work (higher energy), there must be a minimum somewhere between
1878 * the old position and the new one.
1880 * Due to the finite numerical accuracy, it turns out that it is a good idea
1881 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1882 * This leads to lower final energies in the tests I've done. / Erik
1887 c = a + stepsize; /* reference position along line is zero */
1889 /* Check stepsize first. We do not allow displacements
1890 * larger than emstep.
1896 for (i = 0; i < n; i++)
1899 if (delta > maxdelta)
1904 if (maxdelta > inputrec->em_stepsize)
1909 while (maxdelta > inputrec->em_stepsize);
1911 /* Take a trial step */
1912 for (i = 0; i < n; i++)
1914 xc[i] = lastx[i] + c*s[i];
1918 /* Calculate energy for the trial step */
1919 ems.s.x = (rvec *)xc;
1921 evaluate_energy(fplog, cr,
1922 top_global, &ems, top,
1923 inputrec, nrnb, wcycle, gstat,
1924 vsite, constr, fcd, graph, mdatoms, fr,
1925 mu_tot, enerd, vir, pres, step, FALSE);
1928 /* Calc derivative along line */
1929 for (gpc = 0, i = 0; i < n; i++)
1931 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1933 /* Sum the gradient along the line across CPUs */
1936 gmx_sumd(1, &gpc, cr);
1939 /* This is the max amount of increase in energy we tolerate */
1940 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1942 /* Accept the step if the energy is lower, or if it is not significantly higher
1943 * and the line derivative is still negative.
1945 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1948 /* Great, we found a better energy. Increase step for next iteration
1949 * if we are still going down, decrease it otherwise
1953 stepsize *= 1.618034; /* The golden section */
1957 stepsize *= 0.618034; /* 1/golden section */
1962 /* New energy is the same or higher. We will have to do some work
1963 * to find a smaller value in the interval. Take smaller step next time!
1966 stepsize *= 0.618034;
1969 /* OK, if we didn't find a lower value we will have to locate one now - there must
1970 * be one in the interval [a=0,c].
1971 * The same thing is valid here, though: Don't spend dozens of iterations to find
1972 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1973 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1975 * I also have a safeguard for potentially really patological functions so we never
1976 * take more than 20 steps before we give up ...
1978 * If we already found a lower value we just skip this step and continue to the update.
1987 /* Select a new trial point.
1988 * If the derivatives at points a & c have different sign we interpolate to zero,
1989 * otherwise just do a bisection.
1992 if (gpa < 0 && gpc > 0)
1994 b = a + gpa*(a-c)/(gpc-gpa);
2001 /* safeguard if interpolation close to machine accuracy causes errors:
2002 * never go outside the interval
2004 if (b <= a || b >= c)
2009 /* Take a trial step */
2010 for (i = 0; i < n; i++)
2012 xb[i] = lastx[i] + b*s[i];
2016 /* Calculate energy for the trial step */
2017 ems.s.x = (rvec *)xb;
2019 evaluate_energy(fplog, cr,
2020 top_global, &ems, top,
2021 inputrec, nrnb, wcycle, gstat,
2022 vsite, constr, fcd, graph, mdatoms, fr,
2023 mu_tot, enerd, vir, pres, step, FALSE);
2028 for (gpb = 0, i = 0; i < n; i++)
2030 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2033 /* Sum the gradient along the line across CPUs */
2036 gmx_sumd(1, &gpb, cr);
2039 /* Keep one of the intervals based on the value of the derivative at the new point */
2042 /* Replace c endpoint with b */
2046 /* swap coord pointers b/c */
2056 /* Replace a endpoint with b */
2060 /* swap coord pointers a/b */
2070 * Stop search as soon as we find a value smaller than the endpoints,
2071 * or if the tolerance is below machine precision.
2072 * Never run more than 20 steps, no matter what.
2076 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2078 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2080 /* OK. We couldn't find a significantly lower energy.
2081 * If ncorr==0 this was steepest descent, and then we give up.
2082 * If not, reset memory to restart as steepest descent before quitting.
2094 /* Search in gradient direction */
2095 for (i = 0; i < n; i++)
2097 dx[point][i] = ff[i];
2099 /* Reset stepsize */
2100 stepsize = 1.0/fnorm;
2105 /* Select min energy state of A & C, put the best in xx/ff/Epot
2111 for (i = 0; i < n; i++)
2122 for (i = 0; i < n; i++)
2136 for (i = 0; i < n; i++)
2144 /* Update the memory information, and calculate a new
2145 * approximation of the inverse hessian
2148 /* Have new data in Epot, xx, ff */
2149 if (ncorr < nmaxcorr)
2154 for (i = 0; i < n; i++)
2156 dg[point][i] = lastf[i]-ff[i];
2157 dx[point][i] *= stepsize;
2162 for (i = 0; i < n; i++)
2164 dgdg += dg[point][i]*dg[point][i];
2165 dgdx += dg[point][i]*dx[point][i];
2170 rho[point] = 1.0/dgdx;
2173 if (point >= nmaxcorr)
2179 for (i = 0; i < n; i++)
2186 /* Recursive update. First go back over the memory points */
2187 for (k = 0; k < ncorr; k++)
2196 for (i = 0; i < n; i++)
2198 sq += dx[cp][i]*p[i];
2201 alpha[cp] = rho[cp]*sq;
2203 for (i = 0; i < n; i++)
2205 p[i] -= alpha[cp]*dg[cp][i];
2209 for (i = 0; i < n; i++)
2214 /* And then go forward again */
2215 for (k = 0; k < ncorr; k++)
2218 for (i = 0; i < n; i++)
2220 yr += p[i]*dg[cp][i];
2224 beta = alpha[cp]-beta;
2226 for (i = 0; i < n; i++)
2228 p[i] += beta*dx[cp][i];
2238 for (i = 0; i < n; i++)
2242 dx[point][i] = p[i];
2252 /* Test whether the convergence criterion is met */
2253 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2255 /* Print it if necessary */
2260 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2261 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2263 /* Store the new (lower) energies */
2264 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2265 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2266 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2267 do_log = do_per_step(step, inputrec->nstlog);
2268 do_ene = do_per_step(step, inputrec->nstenergy);
2271 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2273 print_ebin(outf->fp_ene, do_ene, FALSE, FALSE,
2274 do_log ? fplog : NULL, step, step, eprNORMAL,
2275 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2278 /* Stop when the maximum force lies below tolerance.
2279 * If we have reached machine precision, converged is already set to true.
2282 converged = converged || (fmax < inputrec->em_tol);
2284 } /* End of the loop */
2288 step--; /* we never took that last step in this case */
2291 if (fmax > inputrec->em_tol)
2295 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2296 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2301 /* If we printed energy and/or logfile last step (which was the last step)
2302 * we don't have to do it again, but otherwise print the final values.
2304 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2306 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2308 if (!do_ene || !do_log) /* Write final energy file entries */
2310 print_ebin(outf->fp_ene, !do_ene, FALSE, FALSE,
2311 !do_log ? fplog : NULL, step, step, eprNORMAL,
2312 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2315 /* Print some stuff... */
2318 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2322 * For accurate normal mode calculation it is imperative that we
2323 * store the last conformation into the full precision binary trajectory.
2325 * However, we should only do it if we did NOT already write this step
2326 * above (which we did if do_x or do_f was true).
2328 do_x = !do_per_step(step, inputrec->nstxout);
2329 do_f = !do_per_step(step, inputrec->nstfout);
2330 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2331 top_global, inputrec, step,
2336 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2337 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2338 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2339 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2341 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2344 finish_em(cr, outf, walltime_accounting, wcycle);
2346 /* To print the actual number of steps we needed somewhere */
2347 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2350 } /* That's all folks */
2353 double do_steep(FILE *fplog, t_commrec *cr,
2354 int nfile, const t_filenm fnm[],
2355 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2356 int gmx_unused nstglobalcomm,
2357 gmx_vsite_t *vsite, gmx_constr_t constr,
2358 int gmx_unused stepout,
2359 t_inputrec *inputrec,
2360 gmx_mtop_t *top_global, t_fcdata *fcd,
2361 t_state *state_global,
2363 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2364 gmx_edsam_t gmx_unused ed,
2366 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2367 gmx_membed_t gmx_unused membed,
2368 real gmx_unused cpt_period, real gmx_unused max_hours,
2369 const char gmx_unused *deviceOptions,
2370 unsigned long gmx_unused Flags,
2371 gmx_walltime_accounting_t walltime_accounting)
2373 const char *SD = "Steepest Descents";
2374 em_state_t *s_min, *s_try;
2376 gmx_localtop_t *top;
2377 gmx_enerdata_t *enerd;
2379 gmx_global_stat_t gstat;
2381 real stepsize, constepsize;
2385 gmx_bool bDone, bAbort, do_x, do_f;
2390 int steps_accepted = 0;
2394 s_min = init_em_state();
2395 s_try = init_em_state();
2397 /* Init em and store the local state in s_try */
2398 init_em(fplog, SD, cr, inputrec,
2399 state_global, top_global, s_try, &top, &f, &f_global,
2400 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2401 nfile, fnm, &outf, &mdebin);
2403 /* Print to log file */
2404 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2406 /* Set variables for stepsize (in nm). This is the largest
2407 * step that we are going to make in any direction.
2409 ustep = inputrec->em_stepsize;
2412 /* Max number of steps */
2413 nsteps = inputrec->nsteps;
2417 /* Print to the screen */
2418 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2422 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2425 /**** HERE STARTS THE LOOP ****
2426 * count is the counter for the number of steps
2427 * bDone will be TRUE when the minimization has converged
2428 * bAbort will be TRUE when nsteps steps have been performed or when
2429 * the stepsize becomes smaller than is reasonable for machine precision
2434 while (!bDone && !bAbort)
2436 bAbort = (nsteps >= 0) && (count == nsteps);
2438 /* set new coordinates, except for first step */
2441 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2442 s_min, stepsize, s_min->f, s_try,
2443 constr, top, nrnb, wcycle, count);
2446 evaluate_energy(fplog, cr,
2447 top_global, s_try, top,
2448 inputrec, nrnb, wcycle, gstat,
2449 vsite, constr, fcd, graph, mdatoms, fr,
2450 mu_tot, enerd, vir, pres, count, count == 0);
2454 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2459 s_min->epot = s_try->epot + 1;
2462 /* Print it if necessary */
2467 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2468 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2469 (s_try->epot < s_min->epot) ? '\n' : '\r');
2472 if (s_try->epot < s_min->epot)
2474 /* Store the new (lower) energies */
2475 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2476 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2477 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2478 print_ebin(outf->fp_ene, TRUE,
2479 do_per_step(steps_accepted, inputrec->nstdisreout),
2480 do_per_step(steps_accepted, inputrec->nstorireout),
2481 fplog, count, count, eprNORMAL, TRUE,
2482 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2487 /* Now if the new energy is smaller than the previous...
2488 * or if this is the first step!
2489 * or if we did random steps!
2492 if ( (count == 0) || (s_try->epot < s_min->epot) )
2496 /* Test whether the convergence criterion is met... */
2497 bDone = (s_try->fmax < inputrec->em_tol);
2499 /* Copy the arrays for force, positions and energy */
2500 /* The 'Min' array always holds the coords and forces of the minimal
2502 swap_em_state(s_min, s_try);
2508 /* Write to trn, if necessary */
2509 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2510 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2511 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2512 top_global, inputrec, count,
2513 s_min, state_global, f_global);
2517 /* If energy is not smaller make the step smaller... */
2520 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2522 /* Reload the old state */
2523 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2524 s_min, top, mdatoms, fr, vsite, constr,
2529 /* Determine new step */
2530 stepsize = ustep/s_min->fmax;
2532 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2534 if (count == nsteps || ustep < 1e-12)
2536 if (count == nsteps || ustep < 1e-6)
2541 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2542 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2548 } /* End of the loop */
2550 /* Print some shit... */
2553 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2555 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2556 top_global, inputrec, count,
2557 s_min, state_global, f_global);
2559 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2563 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2564 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2565 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2566 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2569 finish_em(cr, outf, walltime_accounting, wcycle);
2571 /* To print the actual number of steps we needed somewhere */
2572 inputrec->nsteps = count;
2574 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2577 } /* That's all folks */
2580 double do_nm(FILE *fplog, t_commrec *cr,
2581 int nfile, const t_filenm fnm[],
2582 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2583 int gmx_unused nstglobalcomm,
2584 gmx_vsite_t *vsite, gmx_constr_t constr,
2585 int gmx_unused stepout,
2586 t_inputrec *inputrec,
2587 gmx_mtop_t *top_global, t_fcdata *fcd,
2588 t_state *state_global,
2590 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2591 gmx_edsam_t gmx_unused ed,
2593 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2594 gmx_membed_t gmx_unused membed,
2595 real gmx_unused cpt_period, real gmx_unused max_hours,
2596 const char gmx_unused *deviceOptions,
2597 unsigned long gmx_unused Flags,
2598 gmx_walltime_accounting_t walltime_accounting)
2600 const char *NM = "Normal Mode Analysis";
2602 int natoms, atom, d;
2605 gmx_localtop_t *top;
2606 gmx_enerdata_t *enerd;
2608 gmx_global_stat_t gstat;
2610 real t, t0, lambda, lam0;
2615 gmx_bool bSparse; /* use sparse matrix storage format */
2617 gmx_sparsematrix_t * sparse_matrix = NULL;
2618 real * full_matrix = NULL;
2619 em_state_t * state_work;
2621 /* added with respect to mdrun */
2622 int i, j, k, row, col;
2623 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2629 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2632 state_work = init_em_state();
2634 /* Init em and store the local state in state_minimum */
2635 init_em(fplog, NM, cr, inputrec,
2636 state_global, top_global, state_work, &top,
2638 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2639 nfile, fnm, &outf, NULL);
2641 natoms = top_global->natoms;
2649 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2650 " which MIGHT not be accurate enough for normal mode analysis.\n"
2651 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2652 " are fairly modest even if you recompile in double precision.\n\n");
2656 /* Check if we can/should use sparse storage format.
2658 * Sparse format is only useful when the Hessian itself is sparse, which it
2659 * will be when we use a cutoff.
2660 * For small systems (n<1000) it is easier to always use full matrix format, though.
2662 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2664 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2667 else if (top_global->natoms < 1000)
2669 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2674 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2680 sz = DIM*top_global->natoms;
2682 fprintf(stderr, "Allocating Hessian memory...\n\n");
2686 sparse_matrix = gmx_sparsematrix_init(sz);
2687 sparse_matrix->compressed_symmetric = TRUE;
2691 snew(full_matrix, sz*sz);
2695 /* Initial values */
2696 t0 = inputrec->init_t;
2697 lam0 = inputrec->fepvals->init_lambda;
2705 /* Write start time and temperature */
2706 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2708 /* fudge nr of steps to nr of atoms */
2709 inputrec->nsteps = natoms*2;
2713 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2714 *(top_global->name), (int)inputrec->nsteps);
2717 nnodes = cr->nnodes;
2719 /* Make evaluate_energy do a single node force calculation */
2721 evaluate_energy(fplog, cr,
2722 top_global, state_work, top,
2723 inputrec, nrnb, wcycle, gstat,
2724 vsite, constr, fcd, graph, mdatoms, fr,
2725 mu_tot, enerd, vir, pres, -1, TRUE);
2726 cr->nnodes = nnodes;
2728 /* if forces are not small, warn user */
2729 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2731 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2732 if (state_work->fmax > 1.0e-3)
2734 md_print_info(cr, fplog,
2735 "The force is probably not small enough to "
2736 "ensure that you are at a minimum.\n"
2737 "Be aware that negative eigenvalues may occur\n"
2738 "when the resulting matrix is diagonalized.\n\n");
2741 /***********************************************************
2743 * Loop over all pairs in matrix
2745 * do_force called twice. Once with positive and
2746 * once with negative displacement
2748 ************************************************************/
2750 /* Steps are divided one by one over the nodes */
2751 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2754 for (d = 0; d < DIM; d++)
2756 x_min = state_work->s.x[atom][d];
2758 state_work->s.x[atom][d] = x_min - der_range;
2760 /* Make evaluate_energy do a single node force calculation */
2762 evaluate_energy(fplog, cr,
2763 top_global, state_work, top,
2764 inputrec, nrnb, wcycle, gstat,
2765 vsite, constr, fcd, graph, mdatoms, fr,
2766 mu_tot, enerd, vir, pres, atom*2, FALSE);
2768 for (i = 0; i < natoms; i++)
2770 copy_rvec(state_work->f[i], fneg[i]);
2773 state_work->s.x[atom][d] = x_min + der_range;
2775 evaluate_energy(fplog, cr,
2776 top_global, state_work, top,
2777 inputrec, nrnb, wcycle, gstat,
2778 vsite, constr, fcd, graph, mdatoms, fr,
2779 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2780 cr->nnodes = nnodes;
2782 /* x is restored to original */
2783 state_work->s.x[atom][d] = x_min;
2785 for (j = 0; j < natoms; j++)
2787 for (k = 0; (k < DIM); k++)
2790 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2798 #define mpi_type MPI_DOUBLE
2800 #define mpi_type MPI_FLOAT
2802 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2803 cr->mpi_comm_mygroup);
2808 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2814 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2815 cr->mpi_comm_mygroup, &stat);
2820 row = (atom + node)*DIM + d;
2822 for (j = 0; j < natoms; j++)
2824 for (k = 0; k < DIM; k++)
2830 if (col >= row && dfdx[j][k] != 0.0)
2832 gmx_sparsematrix_increment_value(sparse_matrix,
2833 row, col, dfdx[j][k]);
2838 full_matrix[row*sz+col] = dfdx[j][k];
2845 if (bVerbose && fplog)
2850 /* write progress */
2851 if (MASTER(cr) && bVerbose)
2853 fprintf(stderr, "\rFinished step %d out of %d",
2854 min(atom+nnodes, natoms), natoms);
2861 fprintf(stderr, "\n\nWriting Hessian...\n");
2862 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2865 finish_em(cr, outf, walltime_accounting, wcycle);
2867 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);