2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file defines integrators for energy minimization
41 * \author Berk Hess <hess@kth.se>
42 * \author Erik Lindahl <erik@kth.se>
43 * \ingroup module_mdrun
56 #include "gromacs/commandline/filenm.h"
57 #include "gromacs/domdec/collect.h"
58 #include "gromacs/domdec/domdec.h"
59 #include "gromacs/domdec/domdec_struct.h"
60 #include "gromacs/domdec/partition.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fileio/confio.h"
63 #include "gromacs/fileio/mtxio.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/imd/imd.h"
67 #include "gromacs/linearalgebra/sparsematrix.h"
68 #include "gromacs/listed-forces/manage-threading.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/vec.h"
71 #include "gromacs/mdlib/constr.h"
72 #include "gromacs/mdlib/force.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/gmx_omp_nthreads.h"
75 #include "gromacs/mdlib/md_support.h"
76 #include "gromacs/mdlib/mdatoms.h"
77 #include "gromacs/mdlib/mdebin.h"
78 #include "gromacs/mdlib/mdrun.h"
79 #include "gromacs/mdlib/mdsetup.h"
80 #include "gromacs/mdlib/ns.h"
81 #include "gromacs/mdlib/shellfc.h"
82 #include "gromacs/mdlib/sim_util.h"
83 #include "gromacs/mdlib/tgroup.h"
84 #include "gromacs/mdlib/trajectory_writing.h"
85 #include "gromacs/mdlib/update.h"
86 #include "gromacs/mdlib/vsite.h"
87 #include "gromacs/mdtypes/commrec.h"
88 #include "gromacs/mdtypes/inputrec.h"
89 #include "gromacs/mdtypes/md_enums.h"
90 #include "gromacs/mdtypes/state.h"
91 #include "gromacs/pbcutil/mshift.h"
92 #include "gromacs/pbcutil/pbc.h"
93 #include "gromacs/timing/wallcycle.h"
94 #include "gromacs/timing/walltime_accounting.h"
95 #include "gromacs/topology/mtop_util.h"
96 #include "gromacs/topology/topology.h"
97 #include "gromacs/utility/cstringutil.h"
98 #include "gromacs/utility/exceptions.h"
99 #include "gromacs/utility/fatalerror.h"
100 #include "gromacs/utility/logger.h"
101 #include "gromacs/utility/smalloc.h"
103 #include "integrator.h"
105 //! Utility structure for manipulating states during EM
107 //! Copy of the global state
110 PaddedVector<gmx::RVec> f;
113 //! Norm of the force
121 //! Print the EM starting conditions
122 static void print_em_start(FILE *fplog,
124 gmx_walltime_accounting_t walltime_accounting,
125 gmx_wallcycle_t wcycle,
128 walltime_accounting_start_time(walltime_accounting);
129 wallcycle_start(wcycle, ewcRUN);
130 print_start(fplog, cr, walltime_accounting, name);
133 //! Stop counting time for EM
134 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
135 gmx_wallcycle_t wcycle)
137 wallcycle_stop(wcycle, ewcRUN);
139 walltime_accounting_end_time(walltime_accounting);
142 //! Printing a log file and console header
143 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
146 fprintf(out, "%s:\n", minimizer);
147 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
148 fprintf(out, " Number of steps = %12d\n", nsteps);
151 //! Print warning message
152 static void warn_step(FILE *fp,
158 constexpr bool realIsDouble = GMX_DOUBLE;
161 if (!std::isfinite(fmax))
164 "\nEnergy minimization has stopped because the force "
165 "on at least one atom is not finite. This usually means "
166 "atoms are overlapping. Modify the input coordinates to "
167 "remove atom overlap or use soft-core potentials with "
168 "the free energy code to avoid infinite forces.\n%s",
170 "You could also be lucky that switching to double precision "
171 "is sufficient to obtain finite forces.\n" :
177 "\nEnergy minimization reached the maximum number "
178 "of steps before the forces reached the requested "
179 "precision Fmax < %g.\n", ftol);
184 "\nEnergy minimization has stopped, but the forces have "
185 "not converged to the requested precision Fmax < %g (which "
186 "may not be possible for your system). It stopped "
187 "because the algorithm tried to make a new step whose size "
188 "was too small, or there was no change in the energy since "
189 "last step. Either way, we regard the minimization as "
190 "converged to within the available machine precision, "
191 "given your starting configuration and EM parameters.\n%s%s",
194 "\nDouble precision normally gives you higher accuracy, but "
195 "this is often not needed for preparing to run molecular "
199 "You might need to increase your constraint accuracy, or turn\n"
200 "off constraints altogether (set constraints = none in mdp file)\n" :
204 fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
205 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
208 //! Print message about convergence of the EM
209 static void print_converged(FILE *fp, const char *alg, real ftol,
210 int64_t count, gmx_bool bDone, int64_t nsteps,
211 const em_state_t *ems, double sqrtNumAtoms)
213 char buf[STEPSTRSIZE];
217 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
218 alg, ftol, gmx_step_str(count, buf));
220 else if (count < nsteps)
222 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
223 "but did not reach the requested Fmax < %g.\n",
224 alg, gmx_step_str(count, buf), ftol);
228 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
229 alg, ftol, gmx_step_str(count, buf));
233 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
234 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
235 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm/sqrtNumAtoms);
237 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
238 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
239 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm/sqrtNumAtoms);
243 //! Compute the norm and max of the force array in parallel
244 static void get_f_norm_max(const t_commrec *cr,
245 t_grpopts *opts, t_mdatoms *mdatoms, const rvec *f,
246 real *fnorm, real *fmax, int *a_fmax)
250 int la_max, a_max, start, end, i, m, gf;
252 /* This routine finds the largest force and returns it.
253 * On parallel machines the global max is taken.
259 end = mdatoms->homenr;
260 if (mdatoms->cFREEZE)
262 for (i = start; i < end; i++)
264 gf = mdatoms->cFREEZE[i];
266 for (m = 0; m < DIM; m++)
268 if (!opts->nFreeze[gf][m])
270 fam += gmx::square(f[i][m]);
283 for (i = start; i < end; i++)
295 if (la_max >= 0 && DOMAINDECOMP(cr))
297 a_max = cr->dd->globalAtomIndices[la_max];
305 snew(sum, 2*cr->nnodes+1);
306 sum[2*cr->nodeid] = fmax2;
307 sum[2*cr->nodeid+1] = a_max;
308 sum[2*cr->nnodes] = fnorm2;
309 gmx_sumd(2*cr->nnodes+1, sum, cr);
310 fnorm2 = sum[2*cr->nnodes];
311 /* Determine the global maximum */
312 for (i = 0; i < cr->nnodes; i++)
314 if (sum[2*i] > fmax2)
317 a_max = gmx::roundToInt(sum[2*i+1]);
325 *fnorm = sqrt(fnorm2);
337 //! Compute the norm of the force
338 static void get_state_f_norm_max(const t_commrec *cr,
339 t_grpopts *opts, t_mdatoms *mdatoms,
342 get_f_norm_max(cr, opts, mdatoms, ems->f.rvec_array(),
343 &ems->fnorm, &ems->fmax, &ems->a_fmax);
346 //! Initialize the energy minimization
347 static void init_em(FILE *fplog,
348 const gmx::MDLogger &mdlog,
351 const gmx_multisim_t *ms,
352 gmx::IMDOutputProvider *outputProvider,
354 const MdrunOptions &mdrunOptions,
355 t_state *state_global, gmx_mtop_t *top_global,
356 em_state_t *ems, gmx_localtop_t **top,
357 t_nrnb *nrnb, rvec mu_tot,
358 t_forcerec *fr, gmx_enerdata_t **enerd,
359 t_graph **graph, gmx::MDAtoms *mdAtoms, gmx_global_stat_t *gstat,
360 gmx_vsite_t *vsite, gmx::Constraints *constr, gmx_shellfc_t **shellfc,
361 int nfile, const t_filenm fnm[],
362 gmx_mdoutf_t *outf, t_mdebin **mdebin,
363 gmx_wallcycle_t wcycle)
369 fprintf(fplog, "Initiating %s\n", title);
374 state_global->ngtc = 0;
376 /* Initialize lambda variables */
377 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, nullptr);
382 /* Interactive molecular dynamics */
383 init_IMD(ir, cr, ms, top_global, fplog, 1,
384 MASTER(cr) ? state_global->x.rvec_array() : nullptr,
385 nfile, fnm, nullptr, mdrunOptions);
389 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
391 *shellfc = init_shell_flexcon(stdout,
393 constr ? constr->numFlexibleConstraints() : 0,
399 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
401 /* With energy minimization, shells and flexible constraints are
402 * automatically minimized when treated like normal DOFS.
404 if (shellfc != nullptr)
410 auto mdatoms = mdAtoms->mdatoms();
411 if (DOMAINDECOMP(cr))
413 *top = dd_init_local_top(top_global);
415 dd_init_local_state(cr->dd, state_global, &ems->s);
417 /* Distribute the charge groups over the nodes from the master node */
418 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
419 state_global, top_global, ir,
420 &ems->s, &ems->f, mdAtoms, *top,
422 nrnb, nullptr, FALSE);
423 dd_store_state(cr->dd, &ems->s);
429 state_change_natoms(state_global, state_global->natoms);
430 /* Just copy the state */
431 ems->s = *state_global;
432 state_change_natoms(&ems->s, ems->s.natoms);
433 ems->f.resizeWithPadding(ems->s.natoms);
436 mdAlgorithmsSetupAtomData(cr, ir, top_global, *top, fr,
438 constr, vsite, shellfc ? *shellfc : nullptr);
442 set_vsite_top(vsite, *top, mdatoms);
446 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
450 // TODO how should this cross-module support dependency be managed?
451 if (ir->eConstrAlg == econtSHAKE &&
452 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
454 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
455 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
458 if (!ir->bContinuation)
460 /* Constrain the starting coordinates */
462 constr->apply(TRUE, TRUE,
464 ems->s.x.rvec_array(),
465 ems->s.x.rvec_array(),
468 ems->s.lambda[efptFEP], &dvdl_constr,
469 nullptr, nullptr, gmx::ConstraintVariable::Positions);
475 *gstat = global_stat_init(ir);
482 *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, nullptr, wcycle);
485 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
488 if (mdebin != nullptr)
490 /* Init bin for energy stuff */
491 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, nullptr);
495 calc_shifts(ems->s.box, fr->shift_vec);
498 //! Finalize the minimization
499 static void finish_em(const t_commrec *cr, gmx_mdoutf_t outf,
500 gmx_walltime_accounting_t walltime_accounting,
501 gmx_wallcycle_t wcycle)
503 if (!thisRankHasDuty(cr, DUTY_PME))
505 /* Tell the PME only node to finish */
506 gmx_pme_send_finish(cr);
511 em_time_end(walltime_accounting, wcycle);
514 //! Swap two different EM states during minimization
515 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
524 //! Save the EM trajectory
525 static void write_em_traj(FILE *fplog, const t_commrec *cr,
527 gmx_bool bX, gmx_bool bF, const char *confout,
528 gmx_mtop_t *top_global,
529 t_inputrec *ir, int64_t step,
531 t_state *state_global,
532 ObservablesHistory *observablesHistory)
538 mdof_flags |= MDOF_X;
542 mdof_flags |= MDOF_F;
545 /* If we want IMD output, set appropriate MDOF flag */
548 mdof_flags |= MDOF_IMD;
551 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
552 top_global, step, static_cast<double>(step),
553 &state->s, state_global, observablesHistory,
556 if (confout != nullptr)
558 if (DOMAINDECOMP(cr))
560 /* If bX=true, x was collected to state_global in the call above */
563 gmx::ArrayRef<gmx::RVec> globalXRef = MASTER(cr) ? makeArrayRef(state_global->x) : gmx::EmptyArrayRef();
564 dd_collect_vec(cr->dd, &state->s, makeArrayRef(state->s.x), globalXRef);
569 /* Copy the local state pointer */
570 state_global = &state->s;
575 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
577 /* Make molecules whole only for confout writing */
578 do_pbc_mtop(fplog, ir->ePBC, state->s.box, top_global,
579 state_global->x.rvec_array());
582 write_sto_conf_mtop(confout,
583 *top_global->name, top_global,
584 state_global->x.rvec_array(), nullptr, ir->ePBC, state->s.box);
589 //! \brief Do one minimization step
591 // \returns true when the step succeeded, false when a constraint error occurred
592 static bool do_em_step(const t_commrec *cr,
593 t_inputrec *ir, t_mdatoms *md,
594 em_state_t *ems1, real a, const PaddedVector<gmx::RVec> *force,
596 gmx::Constraints *constr,
603 int nthreads gmx_unused;
605 bool validStep = true;
610 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
612 gmx_incons("state mismatch in do_em_step");
615 s2->flags = s1->flags;
617 if (s2->natoms != s1->natoms)
619 state_change_natoms(s2, s1->natoms);
620 ems2->f.resizeWithPadding(s2->natoms);
622 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
624 s2->cg_gl.resize(s1->cg_gl.size());
627 copy_mat(s1->box, s2->box);
628 /* Copy free energy state */
629 s2->lambda = s1->lambda;
630 copy_mat(s1->box, s2->box);
635 nthreads = gmx_omp_nthreads_get(emntUpdate);
636 #pragma omp parallel num_threads(nthreads)
638 const rvec *x1 = s1->x.rvec_array();
639 rvec *x2 = s2->x.rvec_array();
640 const rvec *f = force->rvec_array();
643 #pragma omp for schedule(static) nowait
644 for (int i = start; i < end; i++)
652 for (int m = 0; m < DIM; m++)
654 if (ir->opts.nFreeze[gf][m])
660 x2[i][m] = x1[i][m] + a*f[i][m];
664 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
667 if (s2->flags & (1<<estCGP))
669 /* Copy the CG p vector */
670 const rvec *p1 = s1->cg_p.rvec_array();
671 rvec *p2 = s2->cg_p.rvec_array();
672 #pragma omp for schedule(static) nowait
673 for (int i = start; i < end; i++)
675 // Trivial OpenMP block that does not throw
676 copy_rvec(p1[i], p2[i]);
680 if (DOMAINDECOMP(cr))
682 /* OpenMP does not supported unsigned loop variables */
683 #pragma omp for schedule(static) nowait
684 for (int i = 0; i < static_cast<int>(s2->cg_gl.size()); i++)
686 s2->cg_gl[i] = s1->cg_gl[i];
691 if (DOMAINDECOMP(cr))
693 s2->ddp_count = s1->ddp_count;
694 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
701 constr->apply(TRUE, TRUE,
703 s1->x.rvec_array(), s2->x.rvec_array(),
705 s2->lambda[efptBONDED], &dvdl_constr,
706 nullptr, nullptr, gmx::ConstraintVariable::Positions);
710 /* This global reduction will affect performance at high
711 * parallelization, but we can not really avoid it.
712 * But usually EM is not run at high parallelization.
714 int reductionBuffer = static_cast<int>(!validStep);
715 gmx_sumi(1, &reductionBuffer, cr);
716 validStep = (reductionBuffer == 0);
719 // We should move this check to the different minimizers
720 if (!validStep && ir->eI != eiSteep)
722 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
723 EI(ir->eI), EI(eiSteep), EI(ir->eI));
730 //! Prepare EM for using domain decomposition parallellization
731 static void em_dd_partition_system(FILE *fplog,
732 const gmx::MDLogger &mdlog,
733 int step, const t_commrec *cr,
734 gmx_mtop_t *top_global, t_inputrec *ir,
735 em_state_t *ems, gmx_localtop_t *top,
736 gmx::MDAtoms *mdAtoms, t_forcerec *fr,
737 gmx_vsite_t *vsite, gmx::Constraints *constr,
738 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
740 /* Repartition the domain decomposition */
741 dd_partition_system(fplog, mdlog, step, cr, FALSE, 1,
742 nullptr, top_global, ir,
744 mdAtoms, top, fr, vsite, constr,
745 nrnb, wcycle, FALSE);
746 dd_store_state(cr->dd, &ems->s);
752 /*! \brief Class to handle the work of setting and doing an energy evaluation.
754 * This class is a mere aggregate of parameters to pass to evaluate an
755 * energy, so that future changes to names and types of them consume
756 * less time when refactoring other code.
758 * Aggregate initialization is used, for which the chief risk is that
759 * if a member is added at the end and not all initializer lists are
760 * updated, then the member will be value initialized, which will
761 * typically mean initialization to zero.
763 * We only want to construct one of these with an initializer list, so
764 * we explicitly delete the default constructor. */
765 class EnergyEvaluator
768 //! We only intend to construct such objects with an initializer list.
769 #if __GNUC__ > 4 || (__GNUC__ == 4 && __GNUC_MINOR__ >= 9)
770 // Aspects of the C++11 spec changed after GCC 4.8.5, and
771 // compilation of the initializer list construction in
772 // runner.cpp fails in GCC 4.8.5.
773 EnergyEvaluator() = delete;
775 /*! \brief Evaluates an energy on the state in \c ems.
777 * \todo In practice, the same objects mu_tot, vir, and pres
778 * are always passed to this function, so we would rather have
779 * them as data members. However, their C-array types are
780 * unsuited for aggregate initialization. When the types
781 * improve, the call signature of this method can be reduced.
783 void run(em_state_t *ems, rvec mu_tot,
784 tensor vir, tensor pres,
785 int64_t count, gmx_bool bFirst);
786 //! Handles logging (deprecated).
789 const gmx::MDLogger &mdlog;
790 //! Handles communication.
792 //! Coordinates multi-simulations.
793 const gmx_multisim_t *ms;
794 //! Holds the simulation topology.
795 gmx_mtop_t *top_global;
796 //! Holds the domain topology.
798 //! User input options.
799 t_inputrec *inputrec;
800 //! Manages flop accounting.
802 //! Manages wall cycle accounting.
803 gmx_wallcycle_t wcycle;
804 //! Coordinates global reduction.
805 gmx_global_stat_t gstat;
806 //! Handles virtual sites.
808 //! Handles constraints.
809 gmx::Constraints *constr;
810 //! Handles strange things.
812 //! Molecular graph for SHAKE.
814 //! Per-atom data for this domain.
815 gmx::MDAtoms *mdAtoms;
816 //! Handles how to calculate the forces.
818 //! Schedule of force-calculation work each step for this task.
819 gmx::PpForceWorkload *ppForceWorkload;
820 //! Stores the computed energies.
821 gmx_enerdata_t *enerd;
825 EnergyEvaluator::run(em_state_t *ems, rvec mu_tot,
826 tensor vir, tensor pres,
827 int64_t count, gmx_bool bFirst)
831 tensor force_vir, shake_vir, ekin;
832 real dvdl_constr, prescorr, enercorr, dvdlcorr;
835 /* Set the time to the initial time, the time does not change during EM */
836 t = inputrec->init_t;
839 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
841 /* This is the first state or an old state used before the last ns */
847 if (inputrec->nstlist > 0)
855 construct_vsites(vsite, ems->s.x.rvec_array(), 1, nullptr,
856 top->idef.iparams, top->idef.il,
857 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
860 if (DOMAINDECOMP(cr) && bNS)
862 /* Repartition the domain decomposition */
863 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec,
864 ems, top, mdAtoms, fr, vsite, constr,
868 /* Calc force & energy on new trial position */
869 /* do_force always puts the charge groups in the box and shifts again
870 * We do not unshift, so molecules are always whole in congrad.c
872 do_force(fplog, cr, ms, inputrec, nullptr, nullptr,
873 count, nrnb, wcycle, top, &top_global->groups,
874 ems->s.box, ems->s.x.arrayRefWithPadding(), &ems->s.hist,
875 ems->f.arrayRefWithPadding(), force_vir, mdAtoms->mdatoms(), enerd, fcd,
876 ems->s.lambda, graph, fr, ppForceWorkload, vsite, mu_tot, t, nullptr,
877 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
878 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
879 (bNS ? GMX_FORCE_NS : 0),
881 DdOpenBalanceRegionBeforeForceComputation::yes :
882 DdOpenBalanceRegionBeforeForceComputation::no,
884 DdCloseBalanceRegionAfterForceComputation::yes :
885 DdCloseBalanceRegionAfterForceComputation::no);
887 /* Clear the unused shake virial and pressure */
888 clear_mat(shake_vir);
891 /* Communicate stuff when parallel */
892 if (PAR(cr) && inputrec->eI != eiNM)
894 wallcycle_start(wcycle, ewcMoveE);
896 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
897 inputrec, nullptr, nullptr, nullptr, 1, &terminate,
903 wallcycle_stop(wcycle, ewcMoveE);
906 /* Calculate long range corrections to pressure and energy */
907 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
908 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
909 enerd->term[F_DISPCORR] = enercorr;
910 enerd->term[F_EPOT] += enercorr;
911 enerd->term[F_PRES] += prescorr;
912 enerd->term[F_DVDL] += dvdlcorr;
914 ems->epot = enerd->term[F_EPOT];
918 /* Project out the constraint components of the force */
920 rvec *f_rvec = ems->f.rvec_array();
921 constr->apply(FALSE, FALSE,
923 ems->s.x.rvec_array(), f_rvec, f_rvec,
925 ems->s.lambda[efptBONDED], &dvdl_constr,
926 nullptr, &shake_vir, gmx::ConstraintVariable::ForceDispl);
927 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
928 m_add(force_vir, shake_vir, vir);
932 copy_mat(force_vir, vir);
936 enerd->term[F_PRES] =
937 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
939 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
941 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
943 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
949 //! Parallel utility summing energies and forces
950 static double reorder_partsum(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
951 gmx_mtop_t *top_global,
952 em_state_t *s_min, em_state_t *s_b)
955 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
957 unsigned char *grpnrFREEZE;
961 fprintf(debug, "Doing reorder_partsum\n");
964 const rvec *fm = s_min->f.rvec_array();
965 const rvec *fb = s_b->f.rvec_array();
967 cgs_gl = dd_charge_groups_global(cr->dd);
968 index = cgs_gl->index;
970 /* Collect fm in a global vector fmg.
971 * This conflicts with the spirit of domain decomposition,
972 * but to fully optimize this a much more complicated algorithm is required.
975 snew(fmg, top_global->natoms);
977 ncg = s_min->s.cg_gl.size();
978 cg_gl = s_min->s.cg_gl.data();
980 for (c = 0; c < ncg; c++)
985 for (a = a0; a < a1; a++)
987 copy_rvec(fm[i], fmg[a]);
991 gmx_sum(top_global->natoms*3, fmg[0], cr);
993 /* Now we will determine the part of the sum for the cgs in state s_b */
994 ncg = s_b->s.cg_gl.size();
995 cg_gl = s_b->s.cg_gl.data();
999 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
1000 for (c = 0; c < ncg; c++)
1005 for (a = a0; a < a1; a++)
1007 if (mdatoms->cFREEZE && grpnrFREEZE)
1009 gf = grpnrFREEZE[i];
1011 for (m = 0; m < DIM; m++)
1013 if (!opts->nFreeze[gf][m])
1015 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
1027 //! Print some stuff, like beta, whatever that means.
1028 static real pr_beta(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
1029 gmx_mtop_t *top_global,
1030 em_state_t *s_min, em_state_t *s_b)
1034 /* This is just the classical Polak-Ribiere calculation of beta;
1035 * it looks a bit complicated since we take freeze groups into account,
1036 * and might have to sum it in parallel runs.
1039 if (!DOMAINDECOMP(cr) ||
1040 (s_min->s.ddp_count == cr->dd->ddp_count &&
1041 s_b->s.ddp_count == cr->dd->ddp_count))
1043 const rvec *fm = s_min->f.rvec_array();
1044 const rvec *fb = s_b->f.rvec_array();
1047 /* This part of code can be incorrect with DD,
1048 * since the atom ordering in s_b and s_min might differ.
1050 for (int i = 0; i < mdatoms->homenr; i++)
1052 if (mdatoms->cFREEZE)
1054 gf = mdatoms->cFREEZE[i];
1056 for (int m = 0; m < DIM; m++)
1058 if (!opts->nFreeze[gf][m])
1060 sum += (fb[i][m] - fm[i][m])*fb[i][m];
1067 /* We need to reorder cgs while summing */
1068 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
1072 gmx_sumd(1, &sum, cr);
1075 return sum/gmx::square(s_min->fnorm);
1084 const char *CG = "Polak-Ribiere Conjugate Gradients";
1086 gmx_localtop_t *top;
1087 gmx_enerdata_t *enerd;
1088 gmx_global_stat_t gstat;
1090 double tmp, minstep;
1092 real a, b, c, beta = 0.0;
1096 gmx_bool converged, foundlower;
1098 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1100 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1102 int m, step, nminstep;
1103 auto mdatoms = mdAtoms->mdatoms();
1105 GMX_LOG(mdlog.info).asParagraph().
1106 appendText("Note that activating conjugate gradient energy minimization via the "
1107 "integrator .mdp option and the command gmx mdrun may "
1108 "be available in a different form in a future version of GROMACS, "
1109 "e.g. gmx minimize and an .mdp option.");
1115 // In CG, the state is extended with a search direction
1116 state_global->flags |= (1<<estCGP);
1118 // Ensure the extra per-atom state array gets allocated
1119 state_change_natoms(state_global, state_global->natoms);
1121 // Initialize the search direction to zero
1122 for (RVec &cg_p : state_global->cg_p)
1128 /* Create 4 states on the stack and extract pointers that we will swap */
1129 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1130 em_state_t *s_min = &s0;
1131 em_state_t *s_a = &s1;
1132 em_state_t *s_b = &s2;
1133 em_state_t *s_c = &s3;
1135 /* Init em and store the local state in s_min */
1136 init_em(fplog, mdlog, CG, cr, ms, outputProvider, inputrec, mdrunOptions,
1137 state_global, top_global, s_min, &top,
1138 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1139 vsite, constr, nullptr,
1140 nfile, fnm, &outf, &mdebin, wcycle);
1142 /* Print to log file */
1143 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1145 /* Max number of steps */
1146 number_steps = inputrec->nsteps;
1150 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1154 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1157 EnergyEvaluator energyEvaluator {
1158 fplog, mdlog, cr, ms,
1160 inputrec, nrnb, wcycle, gstat,
1161 vsite, constr, fcd, graph,
1162 mdAtoms, fr, ppForceWorkload, enerd
1164 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1165 /* do_force always puts the charge groups in the box and shifts again
1166 * We do not unshift, so molecules are always whole in congrad.c
1168 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1172 /* Copy stuff to the energy bin for easy printing etc. */
1173 matrix nullBox = {};
1174 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(step),
1175 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1176 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1178 print_ebin_header(fplog, step, step);
1179 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1180 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1183 /* Estimate/guess the initial stepsize */
1184 stepsize = inputrec->em_stepsize/s_min->fnorm;
1188 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1189 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1190 s_min->fmax, s_min->a_fmax+1);
1191 fprintf(stderr, " F-Norm = %12.5e\n",
1192 s_min->fnorm/sqrtNumAtoms);
1193 fprintf(stderr, "\n");
1194 /* and copy to the log file too... */
1195 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1196 s_min->fmax, s_min->a_fmax+1);
1197 fprintf(fplog, " F-Norm = %12.5e\n",
1198 s_min->fnorm/sqrtNumAtoms);
1199 fprintf(fplog, "\n");
1201 /* Start the loop over CG steps.
1202 * Each successful step is counted, and we continue until
1203 * we either converge or reach the max number of steps.
1206 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1209 /* start taking steps in a new direction
1210 * First time we enter the routine, beta=0, and the direction is
1211 * simply the negative gradient.
1214 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1215 rvec *pm = s_min->s.cg_p.rvec_array();
1216 const rvec *sfm = s_min->f.rvec_array();
1219 for (int i = 0; i < mdatoms->homenr; i++)
1221 if (mdatoms->cFREEZE)
1223 gf = mdatoms->cFREEZE[i];
1225 for (m = 0; m < DIM; m++)
1227 if (!inputrec->opts.nFreeze[gf][m])
1229 pm[i][m] = sfm[i][m] + beta*pm[i][m];
1230 gpa -= pm[i][m]*sfm[i][m];
1231 /* f is negative gradient, thus the sign */
1240 /* Sum the gradient along the line across CPUs */
1243 gmx_sumd(1, &gpa, cr);
1246 /* Calculate the norm of the search vector */
1247 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1249 /* Just in case stepsize reaches zero due to numerical precision... */
1252 stepsize = inputrec->em_stepsize/pnorm;
1256 * Double check the value of the derivative in the search direction.
1257 * If it is positive it must be due to the old information in the
1258 * CG formula, so just remove that and start over with beta=0.
1259 * This corresponds to a steepest descent step.
1264 step--; /* Don't count this step since we are restarting */
1265 continue; /* Go back to the beginning of the big for-loop */
1268 /* Calculate minimum allowed stepsize, before the average (norm)
1269 * relative change in coordinate is smaller than precision
1272 auto s_min_x = makeArrayRef(s_min->s.x);
1273 for (int i = 0; i < mdatoms->homenr; i++)
1275 for (m = 0; m < DIM; m++)
1277 tmp = fabs(s_min_x[i][m]);
1286 /* Add up from all CPUs */
1289 gmx_sumd(1, &minstep, cr);
1292 minstep = GMX_REAL_EPS/sqrt(minstep/(3*top_global->natoms));
1294 if (stepsize < minstep)
1300 /* Write coordinates if necessary */
1301 do_x = do_per_step(step, inputrec->nstxout);
1302 do_f = do_per_step(step, inputrec->nstfout);
1304 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1305 top_global, inputrec, step,
1306 s_min, state_global, observablesHistory);
1308 /* Take a step downhill.
1309 * In theory, we should minimize the function along this direction.
1310 * That is quite possible, but it turns out to take 5-10 function evaluations
1311 * for each line. However, we dont really need to find the exact minimum -
1312 * it is much better to start a new CG step in a modified direction as soon
1313 * as we are close to it. This will save a lot of energy evaluations.
1315 * In practice, we just try to take a single step.
1316 * If it worked (i.e. lowered the energy), we increase the stepsize but
1317 * the continue straight to the next CG step without trying to find any minimum.
1318 * If it didn't work (higher energy), there must be a minimum somewhere between
1319 * the old position and the new one.
1321 * Due to the finite numerical accuracy, it turns out that it is a good idea
1322 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1323 * This leads to lower final energies in the tests I've done. / Erik
1325 s_a->epot = s_min->epot;
1327 c = a + stepsize; /* reference position along line is zero */
1329 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1331 em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec,
1332 s_min, top, mdAtoms, fr, vsite, constr,
1336 /* Take a trial step (new coords in s_c) */
1337 do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c,
1341 /* Calculate energy for the trial step */
1342 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1344 /* Calc derivative along line */
1345 const rvec *pc = s_c->s.cg_p.rvec_array();
1346 const rvec *sfc = s_c->f.rvec_array();
1348 for (int i = 0; i < mdatoms->homenr; i++)
1350 for (m = 0; m < DIM; m++)
1352 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1355 /* Sum the gradient along the line across CPUs */
1358 gmx_sumd(1, &gpc, cr);
1361 /* This is the max amount of increase in energy we tolerate */
1362 tmp = std::sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1364 /* Accept the step if the energy is lower, or if it is not significantly higher
1365 * and the line derivative is still negative.
1367 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1370 /* Great, we found a better energy. Increase step for next iteration
1371 * if we are still going down, decrease it otherwise
1375 stepsize *= 1.618034; /* The golden section */
1379 stepsize *= 0.618034; /* 1/golden section */
1384 /* New energy is the same or higher. We will have to do some work
1385 * to find a smaller value in the interval. Take smaller step next time!
1388 stepsize *= 0.618034;
1394 /* OK, if we didn't find a lower value we will have to locate one now - there must
1395 * be one in the interval [a=0,c].
1396 * The same thing is valid here, though: Don't spend dozens of iterations to find
1397 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1398 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1400 * I also have a safeguard for potentially really pathological functions so we never
1401 * take more than 20 steps before we give up ...
1403 * If we already found a lower value we just skip this step and continue to the update.
1412 /* Select a new trial point.
1413 * If the derivatives at points a & c have different sign we interpolate to zero,
1414 * otherwise just do a bisection.
1416 if (gpa < 0 && gpc > 0)
1418 b = a + gpa*(a-c)/(gpc-gpa);
1425 /* safeguard if interpolation close to machine accuracy causes errors:
1426 * never go outside the interval
1428 if (b <= a || b >= c)
1433 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1435 /* Reload the old state */
1436 em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec,
1437 s_min, top, mdAtoms, fr, vsite, constr,
1441 /* Take a trial step to this new point - new coords in s_b */
1442 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b,
1446 /* Calculate energy for the trial step */
1447 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1449 /* p does not change within a step, but since the domain decomposition
1450 * might change, we have to use cg_p of s_b here.
1452 const rvec *pb = s_b->s.cg_p.rvec_array();
1453 const rvec *sfb = s_b->f.rvec_array();
1455 for (int i = 0; i < mdatoms->homenr; i++)
1457 for (m = 0; m < DIM; m++)
1459 gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1462 /* Sum the gradient along the line across CPUs */
1465 gmx_sumd(1, &gpb, cr);
1470 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1471 s_a->epot, s_b->epot, s_c->epot, gpb);
1474 epot_repl = s_b->epot;
1476 /* Keep one of the intervals based on the value of the derivative at the new point */
1479 /* Replace c endpoint with b */
1480 swap_em_state(&s_b, &s_c);
1486 /* Replace a endpoint with b */
1487 swap_em_state(&s_b, &s_a);
1493 * Stop search as soon as we find a value smaller than the endpoints.
1494 * Never run more than 20 steps, no matter what.
1498 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1501 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1504 /* OK. We couldn't find a significantly lower energy.
1505 * If beta==0 this was steepest descent, and then we give up.
1506 * If not, set beta=0 and restart with steepest descent before quitting.
1516 /* Reset memory before giving up */
1522 /* Select min energy state of A & C, put the best in B.
1524 if (s_c->epot < s_a->epot)
1528 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1529 s_c->epot, s_a->epot);
1531 swap_em_state(&s_b, &s_c);
1538 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1539 s_a->epot, s_c->epot);
1541 swap_em_state(&s_b, &s_a);
1550 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1553 swap_em_state(&s_b, &s_c);
1557 /* new search direction */
1558 /* beta = 0 means forget all memory and restart with steepest descents. */
1559 if (nstcg && ((step % nstcg) == 0))
1565 /* s_min->fnorm cannot be zero, because then we would have converged
1569 /* Polak-Ribiere update.
1570 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1572 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1574 /* Limit beta to prevent oscillations */
1575 if (fabs(beta) > 5.0)
1581 /* update positions */
1582 swap_em_state(&s_min, &s_b);
1585 /* Print it if necessary */
1588 if (mdrunOptions.verbose)
1590 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1591 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1592 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1593 s_min->fmax, s_min->a_fmax+1);
1596 /* Store the new (lower) energies */
1597 matrix nullBox = {};
1598 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(step),
1599 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1600 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1602 do_log = do_per_step(step, inputrec->nstlog);
1603 do_ene = do_per_step(step, inputrec->nstenergy);
1605 /* Prepare IMD energy record, if bIMD is TRUE. */
1606 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1610 print_ebin_header(fplog, step, step);
1612 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1613 do_log ? fplog : nullptr, step, step, eprNORMAL,
1614 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1617 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1618 if (MASTER(cr) && do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x.rvec_array(), inputrec, 0, wcycle))
1620 IMD_send_positions(inputrec->imd);
1623 /* Stop when the maximum force lies below tolerance.
1624 * If we have reached machine precision, converged is already set to true.
1626 converged = converged || (s_min->fmax < inputrec->em_tol);
1628 } /* End of the loop */
1630 /* IMD cleanup, if bIMD is TRUE. */
1631 IMD_finalize(inputrec->bIMD, inputrec->imd);
1635 step--; /* we never took that last step in this case */
1638 if (s_min->fmax > inputrec->em_tol)
1642 warn_step(fplog, inputrec->em_tol, s_min->fmax,
1643 step-1 == number_steps, FALSE);
1650 /* If we printed energy and/or logfile last step (which was the last step)
1651 * we don't have to do it again, but otherwise print the final values.
1655 /* Write final value to log since we didn't do anything the last step */
1656 print_ebin_header(fplog, step, step);
1658 if (!do_ene || !do_log)
1660 /* Write final energy file entries */
1661 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1662 !do_log ? fplog : nullptr, step, step, eprNORMAL,
1663 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1667 /* Print some stuff... */
1670 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1674 * For accurate normal mode calculation it is imperative that we
1675 * store the last conformation into the full precision binary trajectory.
1677 * However, we should only do it if we did NOT already write this step
1678 * above (which we did if do_x or do_f was true).
1680 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1681 * in the trajectory with the same step number.
1683 do_x = !do_per_step(step, inputrec->nstxout);
1684 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1686 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1687 top_global, inputrec, step,
1688 s_min, state_global, observablesHistory);
1693 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1694 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1695 s_min, sqrtNumAtoms);
1696 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1697 s_min, sqrtNumAtoms);
1699 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1702 finish_em(cr, outf, walltime_accounting, wcycle);
1704 /* To print the actual number of steps we needed somewhere */
1705 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1710 Integrator::do_lbfgs()
1712 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1714 gmx_localtop_t *top;
1715 gmx_enerdata_t *enerd;
1716 gmx_global_stat_t gstat;
1718 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1719 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1720 real *rho, *alpha, *p, *s, **dx, **dg;
1721 real a, b, c, maxdelta, delta;
1723 real dgdx, dgdg, sq, yr, beta;
1727 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1729 int start, end, number_steps;
1731 int i, k, m, n, gf, step;
1733 auto mdatoms = mdAtoms->mdatoms();
1735 GMX_LOG(mdlog.info).asParagraph().
1736 appendText("Note that activating L-BFGS energy minimization via the "
1737 "integrator .mdp option and the command gmx mdrun may "
1738 "be available in a different form in a future version of GROMACS, "
1739 "e.g. gmx minimize and an .mdp option.");
1743 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1746 if (nullptr != constr)
1748 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).");
1751 n = 3*state_global->natoms;
1752 nmaxcorr = inputrec->nbfgscorr;
1757 snew(rho, nmaxcorr);
1758 snew(alpha, nmaxcorr);
1761 for (i = 0; i < nmaxcorr; i++)
1767 for (i = 0; i < nmaxcorr; i++)
1776 init_em(fplog, mdlog, LBFGS, cr, ms, outputProvider, inputrec, mdrunOptions,
1777 state_global, top_global, &ems, &top,
1778 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1779 vsite, constr, nullptr,
1780 nfile, fnm, &outf, &mdebin, wcycle);
1783 end = mdatoms->homenr;
1785 /* We need 4 working states */
1786 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1787 em_state_t *sa = &s0;
1788 em_state_t *sb = &s1;
1789 em_state_t *sc = &s2;
1790 em_state_t *last = &s3;
1791 /* Initialize by copying the state from ems (we could skip x and f here) */
1796 /* Print to log file */
1797 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1799 do_log = do_ene = do_x = do_f = TRUE;
1801 /* Max number of steps */
1802 number_steps = inputrec->nsteps;
1804 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1806 for (i = start; i < end; i++)
1808 if (mdatoms->cFREEZE)
1810 gf = mdatoms->cFREEZE[i];
1812 for (m = 0; m < DIM; m++)
1814 frozen[3*i+m] = (inputrec->opts.nFreeze[gf][m] != 0);
1819 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1823 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1828 construct_vsites(vsite, state_global->x.rvec_array(), 1, nullptr,
1829 top->idef.iparams, top->idef.il,
1830 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1833 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1834 /* do_force always puts the charge groups in the box and shifts again
1835 * We do not unshift, so molecules are always whole
1838 EnergyEvaluator energyEvaluator {
1839 fplog, mdlog, cr, ms,
1841 inputrec, nrnb, wcycle, gstat,
1842 vsite, constr, fcd, graph,
1843 mdAtoms, fr, ppForceWorkload, enerd
1845 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1849 /* Copy stuff to the energy bin for easy printing etc. */
1850 matrix nullBox = {};
1851 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(step),
1852 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1853 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1855 print_ebin_header(fplog, step, step);
1856 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1857 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1860 /* Set the initial step.
1861 * since it will be multiplied by the non-normalized search direction
1862 * vector (force vector the first time), we scale it by the
1863 * norm of the force.
1868 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1869 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1870 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1871 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1872 fprintf(stderr, "\n");
1873 /* and copy to the log file too... */
1874 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1875 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1876 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1877 fprintf(fplog, "\n");
1880 // Point is an index to the memory of search directions, where 0 is the first one.
1883 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1884 real *fInit = static_cast<real *>(ems.f.rvec_array()[0]);
1885 for (i = 0; i < n; i++)
1889 dx[point][i] = fInit[i]; /* Initial search direction */
1897 // Stepsize will be modified during the search, and actually it is not critical
1898 // (the main efficiency in the algorithm comes from changing directions), but
1899 // we still need an initial value, so estimate it as the inverse of the norm
1900 // so we take small steps where the potential fluctuates a lot.
1901 stepsize = 1.0/ems.fnorm;
1903 /* Start the loop over BFGS steps.
1904 * Each successful step is counted, and we continue until
1905 * we either converge or reach the max number of steps.
1910 /* Set the gradient from the force */
1912 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1915 /* Write coordinates if necessary */
1916 do_x = do_per_step(step, inputrec->nstxout);
1917 do_f = do_per_step(step, inputrec->nstfout);
1922 mdof_flags |= MDOF_X;
1927 mdof_flags |= MDOF_F;
1932 mdof_flags |= MDOF_IMD;
1935 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1936 top_global, step, static_cast<real>(step), &ems.s, state_global, observablesHistory, ems.f);
1938 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1940 /* make s a pointer to current search direction - point=0 first time we get here */
1943 real *xx = static_cast<real *>(ems.s.x.rvec_array()[0]);
1944 real *ff = static_cast<real *>(ems.f.rvec_array()[0]);
1946 // calculate line gradient in position A
1947 for (gpa = 0, i = 0; i < n; i++)
1952 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1953 * relative change in coordinate is smaller than precision
1955 for (minstep = 0, i = 0; i < n; i++)
1965 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1967 if (stepsize < minstep)
1973 // Before taking any steps along the line, store the old position
1975 real *lastx = static_cast<real *>(last->s.x.data()[0]);
1976 real *lastf = static_cast<real *>(last->f.data()[0]);
1981 /* Take a step downhill.
1982 * In theory, we should find the actual minimum of the function in this
1983 * direction, somewhere along the line.
1984 * That is quite possible, but it turns out to take 5-10 function evaluations
1985 * for each line. However, we dont really need to find the exact minimum -
1986 * it is much better to start a new BFGS step in a modified direction as soon
1987 * as we are close to it. This will save a lot of energy evaluations.
1989 * In practice, we just try to take a single step.
1990 * If it worked (i.e. lowered the energy), we increase the stepsize but
1991 * continue straight to the next BFGS step without trying to find any minimum,
1992 * i.e. we change the search direction too. If the line was smooth, it is
1993 * likely we are in a smooth region, and then it makes sense to take longer
1994 * steps in the modified search direction too.
1996 * If it didn't work (higher energy), there must be a minimum somewhere between
1997 * the old position and the new one. Then we need to start by finding a lower
1998 * value before we change search direction. Since the energy was apparently
1999 * quite rough, we need to decrease the step size.
2001 * Due to the finite numerical accuracy, it turns out that it is a good idea
2002 * to accept a SMALL increase in energy, if the derivative is still downhill.
2003 * This leads to lower final energies in the tests I've done. / Erik
2006 // State "A" is the first position along the line.
2007 // reference position along line is initially zero
2010 // Check stepsize first. We do not allow displacements
2011 // larger than emstep.
2015 // Pick a new position C by adding stepsize to A.
2018 // Calculate what the largest change in any individual coordinate
2019 // would be (translation along line * gradient along line)
2021 for (i = 0; i < n; i++)
2024 if (delta > maxdelta)
2029 // If any displacement is larger than the stepsize limit, reduce the step
2030 if (maxdelta > inputrec->em_stepsize)
2035 while (maxdelta > inputrec->em_stepsize);
2037 // Take a trial step and move the coordinate array xc[] to position C
2038 real *xc = static_cast<real *>(sc->s.x.rvec_array()[0]);
2039 for (i = 0; i < n; i++)
2041 xc[i] = lastx[i] + c*s[i];
2045 // Calculate energy for the trial step in position C
2046 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2048 // Calc line gradient in position C
2049 real *fc = static_cast<real *>(sc->f.rvec_array()[0]);
2050 for (gpc = 0, i = 0; i < n; i++)
2052 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2054 /* Sum the gradient along the line across CPUs */
2057 gmx_sumd(1, &gpc, cr);
2060 // This is the max amount of increase in energy we tolerate.
2061 // By allowing VERY small changes (close to numerical precision) we
2062 // frequently find even better (lower) final energies.
2063 tmp = std::sqrt(GMX_REAL_EPS)*fabs(sa->epot);
2065 // Accept the step if the energy is lower in the new position C (compared to A),
2066 // or if it is not significantly higher and the line derivative is still negative.
2067 foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2068 // If true, great, we found a better energy. We no longer try to alter the
2069 // stepsize, but simply accept this new better position. The we select a new
2070 // search direction instead, which will be much more efficient than continuing
2071 // to take smaller steps along a line. Set fnorm based on the new C position,
2072 // which will be used to update the stepsize to 1/fnorm further down.
2074 // If false, the energy is NOT lower in point C, i.e. it will be the same
2075 // or higher than in point A. In this case it is pointless to move to point C,
2076 // so we will have to do more iterations along the same line to find a smaller
2077 // value in the interval [A=0.0,C].
2078 // Here, A is still 0.0, but that will change when we do a search in the interval
2079 // [0.0,C] below. That search we will do by interpolation or bisection rather
2080 // than with the stepsize, so no need to modify it. For the next search direction
2081 // it will be reset to 1/fnorm anyway.
2085 // OK, if we didn't find a lower value we will have to locate one now - there must
2086 // be one in the interval [a,c].
2087 // The same thing is valid here, though: Don't spend dozens of iterations to find
2088 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2089 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2090 // I also have a safeguard for potentially really pathological functions so we never
2091 // take more than 20 steps before we give up.
2092 // If we already found a lower value we just skip this step and continue to the update.
2097 // Select a new trial point B in the interval [A,C].
2098 // If the derivatives at points a & c have different sign we interpolate to zero,
2099 // otherwise just do a bisection since there might be multiple minima/maxima
2100 // inside the interval.
2101 if (gpa < 0 && gpc > 0)
2103 b = a + gpa*(a-c)/(gpc-gpa);
2110 /* safeguard if interpolation close to machine accuracy causes errors:
2111 * never go outside the interval
2113 if (b <= a || b >= c)
2118 // Take a trial step to point B
2119 real *xb = static_cast<real *>(sb->s.x.rvec_array()[0]);
2120 for (i = 0; i < n; i++)
2122 xb[i] = lastx[i] + b*s[i];
2126 // Calculate energy for the trial step in point B
2127 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2130 // Calculate gradient in point B
2131 real *fb = static_cast<real *>(sb->f.rvec_array()[0]);
2132 for (gpb = 0, i = 0; i < n; i++)
2134 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2137 /* Sum the gradient along the line across CPUs */
2140 gmx_sumd(1, &gpb, cr);
2143 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2144 // at the new point B, and rename the endpoints of this new interval A and C.
2147 /* Replace c endpoint with b */
2149 /* copy state b to c */
2154 /* Replace a endpoint with b */
2156 /* copy state b to a */
2161 * Stop search as soon as we find a value smaller than the endpoints,
2162 * or if the tolerance is below machine precision.
2163 * Never run more than 20 steps, no matter what.
2167 while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2169 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2171 /* OK. We couldn't find a significantly lower energy.
2172 * If ncorr==0 this was steepest descent, and then we give up.
2173 * If not, reset memory to restart as steepest descent before quitting.
2185 /* Search in gradient direction */
2186 for (i = 0; i < n; i++)
2188 dx[point][i] = ff[i];
2190 /* Reset stepsize */
2191 stepsize = 1.0/fnorm;
2196 /* Select min energy state of A & C, put the best in xx/ff/Epot
2198 if (sc->epot < sa->epot)
2220 /* Update the memory information, and calculate a new
2221 * approximation of the inverse hessian
2224 /* Have new data in Epot, xx, ff */
2225 if (ncorr < nmaxcorr)
2230 for (i = 0; i < n; i++)
2232 dg[point][i] = lastf[i]-ff[i];
2233 dx[point][i] *= step_taken;
2238 for (i = 0; i < n; i++)
2240 dgdg += dg[point][i]*dg[point][i];
2241 dgdx += dg[point][i]*dx[point][i];
2246 rho[point] = 1.0/dgdx;
2249 if (point >= nmaxcorr)
2255 for (i = 0; i < n; i++)
2262 /* Recursive update. First go back over the memory points */
2263 for (k = 0; k < ncorr; k++)
2272 for (i = 0; i < n; i++)
2274 sq += dx[cp][i]*p[i];
2277 alpha[cp] = rho[cp]*sq;
2279 for (i = 0; i < n; i++)
2281 p[i] -= alpha[cp]*dg[cp][i];
2285 for (i = 0; i < n; i++)
2290 /* And then go forward again */
2291 for (k = 0; k < ncorr; k++)
2294 for (i = 0; i < n; i++)
2296 yr += p[i]*dg[cp][i];
2300 beta = alpha[cp]-beta;
2302 for (i = 0; i < n; i++)
2304 p[i] += beta*dx[cp][i];
2314 for (i = 0; i < n; i++)
2318 dx[point][i] = p[i];
2326 /* Print it if necessary */
2329 if (mdrunOptions.verbose)
2331 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2332 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2333 step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2336 /* Store the new (lower) energies */
2337 matrix nullBox = {};
2338 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(step),
2339 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2340 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2341 do_log = do_per_step(step, inputrec->nstlog);
2342 do_ene = do_per_step(step, inputrec->nstenergy);
2345 print_ebin_header(fplog, step, step);
2347 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2348 do_log ? fplog : nullptr, step, step, eprNORMAL,
2349 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2352 /* Send x and E to IMD client, if bIMD is TRUE. */
2353 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x.rvec_array(), inputrec, 0, wcycle) && MASTER(cr))
2355 IMD_send_positions(inputrec->imd);
2358 // Reset stepsize in we are doing more iterations
2361 /* Stop when the maximum force lies below tolerance.
2362 * If we have reached machine precision, converged is already set to true.
2364 converged = converged || (ems.fmax < inputrec->em_tol);
2366 } /* End of the loop */
2368 /* IMD cleanup, if bIMD is TRUE. */
2369 IMD_finalize(inputrec->bIMD, inputrec->imd);
2373 step--; /* we never took that last step in this case */
2376 if (ems.fmax > inputrec->em_tol)
2380 warn_step(fplog, inputrec->em_tol, ems.fmax,
2381 step-1 == number_steps, FALSE);
2386 /* If we printed energy and/or logfile last step (which was the last step)
2387 * we don't have to do it again, but otherwise print the final values.
2389 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2391 print_ebin_header(fplog, step, step);
2393 if (!do_ene || !do_log) /* Write final energy file entries */
2395 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2396 !do_log ? fplog : nullptr, step, step, eprNORMAL,
2397 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2400 /* Print some stuff... */
2403 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2407 * For accurate normal mode calculation it is imperative that we
2408 * store the last conformation into the full precision binary trajectory.
2410 * However, we should only do it if we did NOT already write this step
2411 * above (which we did if do_x or do_f was true).
2413 do_x = !do_per_step(step, inputrec->nstxout);
2414 do_f = !do_per_step(step, inputrec->nstfout);
2415 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2416 top_global, inputrec, step,
2417 &ems, state_global, observablesHistory);
2421 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2422 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2423 number_steps, &ems, sqrtNumAtoms);
2424 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2425 number_steps, &ems, sqrtNumAtoms);
2427 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2430 finish_em(cr, outf, walltime_accounting, wcycle);
2432 /* To print the actual number of steps we needed somewhere */
2433 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2437 Integrator::do_steep()
2439 const char *SD = "Steepest Descents";
2440 gmx_localtop_t *top;
2441 gmx_enerdata_t *enerd;
2442 gmx_global_stat_t gstat;
2448 gmx_bool bDone, bAbort, do_x, do_f;
2453 int steps_accepted = 0;
2454 auto mdatoms = mdAtoms->mdatoms();
2456 GMX_LOG(mdlog.info).asParagraph().
2457 appendText("Note that activating steepest-descent energy minimization via the "
2458 "integrator .mdp option and the command gmx mdrun may "
2459 "be available in a different form in a future version of GROMACS, "
2460 "e.g. gmx minimize and an .mdp option.");
2462 /* Create 2 states on the stack and extract pointers that we will swap */
2463 em_state_t s0 {}, s1 {};
2464 em_state_t *s_min = &s0;
2465 em_state_t *s_try = &s1;
2467 /* Init em and store the local state in s_try */
2468 init_em(fplog, mdlog, SD, cr, ms, outputProvider, inputrec, mdrunOptions,
2469 state_global, top_global, s_try, &top,
2470 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2471 vsite, constr, nullptr,
2472 nfile, fnm, &outf, &mdebin, wcycle);
2474 /* Print to log file */
2475 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2477 /* Set variables for stepsize (in nm). This is the largest
2478 * step that we are going to make in any direction.
2480 ustep = inputrec->em_stepsize;
2483 /* Max number of steps */
2484 nsteps = inputrec->nsteps;
2488 /* Print to the screen */
2489 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2493 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2495 EnergyEvaluator energyEvaluator {
2496 fplog, mdlog, cr, ms,
2498 inputrec, nrnb, wcycle, gstat,
2499 vsite, constr, fcd, graph,
2500 mdAtoms, fr, ppForceWorkload, enerd
2503 /**** HERE STARTS THE LOOP ****
2504 * count is the counter for the number of steps
2505 * bDone will be TRUE when the minimization has converged
2506 * bAbort will be TRUE when nsteps steps have been performed or when
2507 * the stepsize becomes smaller than is reasonable for machine precision
2512 while (!bDone && !bAbort)
2514 bAbort = (nsteps >= 0) && (count == nsteps);
2516 /* set new coordinates, except for first step */
2517 bool validStep = true;
2521 do_em_step(cr, inputrec, mdatoms,
2522 s_min, stepsize, &s_min->f, s_try,
2528 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2532 // Signal constraint error during stepping with energy=inf
2533 s_try->epot = std::numeric_limits<real>::infinity();
2538 print_ebin_header(fplog, count, count);
2543 s_min->epot = s_try->epot;
2546 /* Print it if necessary */
2549 if (mdrunOptions.verbose)
2551 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2552 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2553 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2557 if ( (count == 0) || (s_try->epot < s_min->epot) )
2559 /* Store the new (lower) energies */
2560 matrix nullBox = {};
2561 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(count),
2562 mdatoms->tmass, enerd, nullptr, nullptr, nullptr,
2563 nullBox, nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2565 /* Prepare IMD energy record, if bIMD is TRUE. */
2566 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2568 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2569 do_per_step(steps_accepted, inputrec->nstdisreout),
2570 do_per_step(steps_accepted, inputrec->nstorireout),
2571 fplog, count, count, eprNORMAL,
2572 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2577 /* Now if the new energy is smaller than the previous...
2578 * or if this is the first step!
2579 * or if we did random steps!
2582 if ( (count == 0) || (s_try->epot < s_min->epot) )
2586 /* Test whether the convergence criterion is met... */
2587 bDone = (s_try->fmax < inputrec->em_tol);
2589 /* Copy the arrays for force, positions and energy */
2590 /* The 'Min' array always holds the coords and forces of the minimal
2592 swap_em_state(&s_min, &s_try);
2598 /* Write to trn, if necessary */
2599 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2600 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2601 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2602 top_global, inputrec, count,
2603 s_min, state_global, observablesHistory);
2607 /* If energy is not smaller make the step smaller... */
2610 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2612 /* Reload the old state */
2613 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec,
2614 s_min, top, mdAtoms, fr, vsite, constr,
2619 // If the force is very small after finishing minimization,
2620 // we risk dividing by zero when calculating the step size.
2621 // So we check first if the minimization has stopped before
2622 // trying to obtain a new step size.
2625 /* Determine new step */
2626 stepsize = ustep/s_min->fmax;
2629 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2631 if (count == nsteps || ustep < 1e-12)
2633 if (count == nsteps || ustep < 1e-6)
2638 warn_step(fplog, inputrec->em_tol, s_min->fmax,
2639 count == nsteps, constr != nullptr);
2644 /* Send IMD energies and positions, if bIMD is TRUE. */
2645 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box,
2646 MASTER(cr) ? state_global->x.rvec_array() : nullptr,
2647 inputrec, 0, wcycle) &&
2650 IMD_send_positions(inputrec->imd);
2654 } /* End of the loop */
2656 /* IMD cleanup, if bIMD is TRUE. */
2657 IMD_finalize(inputrec->bIMD, inputrec->imd);
2659 /* Print some data... */
2662 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2664 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2665 top_global, inputrec, count,
2666 s_min, state_global, observablesHistory);
2670 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2672 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2673 s_min, sqrtNumAtoms);
2674 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2675 s_min, sqrtNumAtoms);
2678 finish_em(cr, outf, walltime_accounting, wcycle);
2680 /* To print the actual number of steps we needed somewhere */
2681 inputrec->nsteps = count;
2683 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2689 const char *NM = "Normal Mode Analysis";
2692 gmx_localtop_t *top;
2693 gmx_enerdata_t *enerd;
2694 gmx_global_stat_t gstat;
2699 gmx_bool bSparse; /* use sparse matrix storage format */
2701 gmx_sparsematrix_t * sparse_matrix = nullptr;
2702 real * full_matrix = nullptr;
2704 /* added with respect to mdrun */
2706 real der_range = 10.0*std::sqrt(GMX_REAL_EPS);
2708 bool bIsMaster = MASTER(cr);
2709 auto mdatoms = mdAtoms->mdatoms();
2711 GMX_LOG(mdlog.info).asParagraph().
2712 appendText("Note that activating normal-mode analysis via the integrator "
2713 ".mdp option and the command gmx mdrun may "
2714 "be available in a different form in a future version of GROMACS, "
2715 "e.g. gmx normal-modes.");
2717 if (constr != nullptr)
2719 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2722 gmx_shellfc_t *shellfc;
2724 em_state_t state_work {};
2726 /* Init em and store the local state in state_minimum */
2727 init_em(fplog, mdlog, NM, cr, ms, outputProvider, inputrec, mdrunOptions,
2728 state_global, top_global, &state_work, &top,
2729 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2730 vsite, constr, &shellfc,
2731 nfile, fnm, &outf, nullptr, wcycle);
2733 std::vector<int> atom_index = get_atom_index(top_global);
2734 std::vector<gmx::RVec> fneg(atom_index.size(), {0, 0, 0});
2735 snew(dfdx, atom_index.size());
2741 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2742 " which MIGHT not be accurate enough for normal mode analysis.\n"
2743 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2744 " are fairly modest even if you recompile in double precision.\n\n");
2748 /* Check if we can/should use sparse storage format.
2750 * Sparse format is only useful when the Hessian itself is sparse, which it
2751 * will be when we use a cutoff.
2752 * For small systems (n<1000) it is easier to always use full matrix format, though.
2754 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2756 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2759 else if (atom_index.size() < 1000)
2761 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2767 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2771 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2772 sz = DIM*atom_index.size();
2774 fprintf(stderr, "Allocating Hessian memory...\n\n");
2778 sparse_matrix = gmx_sparsematrix_init(sz);
2779 sparse_matrix->compressed_symmetric = TRUE;
2783 snew(full_matrix, sz*sz);
2789 /* Write start time and temperature */
2790 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2792 /* fudge nr of steps to nr of atoms */
2793 inputrec->nsteps = atom_index.size()*2;
2797 fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2798 *(top_global->name), inputrec->nsteps);
2801 nnodes = cr->nnodes;
2803 /* Make evaluate_energy do a single node force calculation */
2805 EnergyEvaluator energyEvaluator {
2806 fplog, mdlog, cr, ms,
2808 inputrec, nrnb, wcycle, gstat,
2809 vsite, constr, fcd, graph,
2810 mdAtoms, fr, ppForceWorkload, enerd
2812 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2813 cr->nnodes = nnodes;
2815 /* if forces are not small, warn user */
2816 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2818 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2819 if (state_work.fmax > 1.0e-3)
2821 GMX_LOG(mdlog.warning).appendText(
2822 "The force is probably not small enough to "
2823 "ensure that you are at a minimum.\n"
2824 "Be aware that negative eigenvalues may occur\n"
2825 "when the resulting matrix is diagonalized.");
2828 /***********************************************************
2830 * Loop over all pairs in matrix
2832 * do_force called twice. Once with positive and
2833 * once with negative displacement
2835 ************************************************************/
2837 /* Steps are divided one by one over the nodes */
2839 auto state_work_x = makeArrayRef(state_work.s.x);
2840 auto state_work_f = makeArrayRef(state_work.f);
2841 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2843 size_t atom = atom_index[aid];
2844 for (size_t d = 0; d < DIM; d++)
2847 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2850 x_min = state_work_x[atom][d];
2852 for (unsigned int dx = 0; (dx < 2); dx++)
2856 state_work_x[atom][d] = x_min - der_range;
2860 state_work_x[atom][d] = x_min + der_range;
2863 /* Make evaluate_energy do a single node force calculation */
2867 /* Now is the time to relax the shells */
2868 relax_shell_flexcon(fplog,
2871 mdrunOptions.verbose,
2882 state_work.f.arrayRefWithPadding(),
2888 &top_global->groups,
2895 DdOpenBalanceRegionBeforeForceComputation::no,
2896 DdCloseBalanceRegionAfterForceComputation::no);
2902 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid*2+dx, FALSE);
2905 cr->nnodes = nnodes;
2909 std::copy(state_work_f.begin(), state_work_f.begin()+atom_index.size(), fneg.begin());
2913 /* x is restored to original */
2914 state_work_x[atom][d] = x_min;
2916 for (size_t j = 0; j < atom_index.size(); j++)
2918 for (size_t k = 0; (k < DIM); k++)
2921 -(state_work_f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2928 #define mpi_type GMX_MPI_REAL
2929 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2930 cr->nodeid, cr->mpi_comm_mygroup);
2935 for (node = 0; (node < nnodes && aid+node < atom_index.size()); node++)
2941 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2942 cr->mpi_comm_mygroup, &stat);
2947 row = (aid + node)*DIM + d;
2949 for (size_t j = 0; j < atom_index.size(); j++)
2951 for (size_t k = 0; k < DIM; k++)
2957 if (col >= row && dfdx[j][k] != 0.0)
2959 gmx_sparsematrix_increment_value(sparse_matrix,
2960 row, col, dfdx[j][k]);
2965 full_matrix[row*sz+col] = dfdx[j][k];
2972 if (mdrunOptions.verbose && fplog)
2977 /* write progress */
2978 if (bIsMaster && mdrunOptions.verbose)
2980 fprintf(stderr, "\rFinished step %d out of %d",
2981 static_cast<int>(std::min(atom+nnodes, atom_index.size())),
2982 static_cast<int>(atom_index.size()));
2989 fprintf(stderr, "\n\nWriting Hessian...\n");
2990 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2993 finish_em(cr, outf, walltime_accounting, wcycle);
2995 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);