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/dlbtiming.h"
59 #include "gromacs/domdec/domdec.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/domdec/mdsetup.h"
62 #include "gromacs/domdec/partition.h"
63 #include "gromacs/ewald/pme.h"
64 #include "gromacs/fileio/confio.h"
65 #include "gromacs/fileio/mtxio.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gmxlib/nrnb.h"
68 #include "gromacs/imd/imd.h"
69 #include "gromacs/linearalgebra/sparsematrix.h"
70 #include "gromacs/listed_forces/manage_threading.h"
71 #include "gromacs/math/functions.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/mdlib/constr.h"
74 #include "gromacs/mdlib/dispersioncorrection.h"
75 #include "gromacs/mdlib/ebin.h"
76 #include "gromacs/mdlib/enerdata_utils.h"
77 #include "gromacs/mdlib/energyoutput.h"
78 #include "gromacs/mdlib/force.h"
79 #include "gromacs/mdlib/forcerec.h"
80 #include "gromacs/mdlib/gmx_omp_nthreads.h"
81 #include "gromacs/mdlib/md_support.h"
82 #include "gromacs/mdlib/mdatoms.h"
83 #include "gromacs/mdlib/stat.h"
84 #include "gromacs/mdlib/tgroup.h"
85 #include "gromacs/mdlib/trajectory_writing.h"
86 #include "gromacs/mdlib/update.h"
87 #include "gromacs/mdlib/vsite.h"
88 #include "gromacs/mdrunutility/handlerestart.h"
89 #include "gromacs/mdrunutility/printtime.h"
90 #include "gromacs/mdtypes/commrec.h"
91 #include "gromacs/mdtypes/inputrec.h"
92 #include "gromacs/mdtypes/md_enums.h"
93 #include "gromacs/mdtypes/mdrunoptions.h"
94 #include "gromacs/mdtypes/state.h"
95 #include "gromacs/pbcutil/mshift.h"
96 #include "gromacs/pbcutil/pbc.h"
97 #include "gromacs/timing/wallcycle.h"
98 #include "gromacs/timing/walltime_accounting.h"
99 #include "gromacs/topology/mtop_util.h"
100 #include "gromacs/topology/topology.h"
101 #include "gromacs/utility/cstringutil.h"
102 #include "gromacs/utility/exceptions.h"
103 #include "gromacs/utility/fatalerror.h"
104 #include "gromacs/utility/logger.h"
105 #include "gromacs/utility/smalloc.h"
107 #include "legacysimulator.h"
110 using gmx::MdrunScheduleWorkload;
112 //! Utility structure for manipulating states during EM
114 //! Copy of the global state
117 PaddedHostVector<gmx::RVec> f;
120 //! Norm of the force
128 //! Print the EM starting conditions
129 static void print_em_start(FILE *fplog,
131 gmx_walltime_accounting_t walltime_accounting,
132 gmx_wallcycle_t wcycle,
135 walltime_accounting_start_time(walltime_accounting);
136 wallcycle_start(wcycle, ewcRUN);
137 print_start(fplog, cr, walltime_accounting, name);
140 //! Stop counting time for EM
141 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
142 gmx_wallcycle_t wcycle)
144 wallcycle_stop(wcycle, ewcRUN);
146 walltime_accounting_end_time(walltime_accounting);
149 //! Printing a log file and console header
150 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
153 fprintf(out, "%s:\n", minimizer);
154 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
155 fprintf(out, " Number of steps = %12d\n", nsteps);
158 //! Print warning message
159 static void warn_step(FILE *fp,
165 constexpr bool realIsDouble = GMX_DOUBLE;
168 if (!std::isfinite(fmax))
171 "\nEnergy minimization has stopped because the force "
172 "on at least one atom is not finite. This usually means "
173 "atoms are overlapping. Modify the input coordinates to "
174 "remove atom overlap or use soft-core potentials with "
175 "the free energy code to avoid infinite forces.\n%s",
177 "You could also be lucky that switching to double precision "
178 "is sufficient to obtain finite forces.\n" :
184 "\nEnergy minimization reached the maximum number "
185 "of steps before the forces reached the requested "
186 "precision Fmax < %g.\n", ftol);
191 "\nEnergy minimization has stopped, but the forces have "
192 "not converged to the requested precision Fmax < %g (which "
193 "may not be possible for your system). It stopped "
194 "because the algorithm tried to make a new step whose size "
195 "was too small, or there was no change in the energy since "
196 "last step. Either way, we regard the minimization as "
197 "converged to within the available machine precision, "
198 "given your starting configuration and EM parameters.\n%s%s",
201 "\nDouble precision normally gives you higher accuracy, but "
202 "this is often not needed for preparing to run molecular "
206 "You might need to increase your constraint accuracy, or turn\n"
207 "off constraints altogether (set constraints = none in mdp file)\n" :
211 fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
212 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
215 //! Print message about convergence of the EM
216 static void print_converged(FILE *fp, const char *alg, real ftol,
217 int64_t count, gmx_bool bDone, int64_t nsteps,
218 const em_state_t *ems, double sqrtNumAtoms)
220 char buf[STEPSTRSIZE];
224 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
225 alg, ftol, gmx_step_str(count, buf));
227 else if (count < nsteps)
229 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
230 "but did not reach the requested Fmax < %g.\n",
231 alg, gmx_step_str(count, buf), ftol);
235 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
236 alg, ftol, gmx_step_str(count, buf));
240 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
241 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
242 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm/sqrtNumAtoms);
244 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
245 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
246 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm/sqrtNumAtoms);
250 //! Compute the norm and max of the force array in parallel
251 static void get_f_norm_max(const t_commrec *cr,
252 t_grpopts *opts, t_mdatoms *mdatoms, const rvec *f,
253 real *fnorm, real *fmax, int *a_fmax)
257 int la_max, a_max, start, end, i, m, gf;
259 /* This routine finds the largest force and returns it.
260 * On parallel machines the global max is taken.
266 end = mdatoms->homenr;
267 if (mdatoms->cFREEZE)
269 for (i = start; i < end; i++)
271 gf = mdatoms->cFREEZE[i];
273 for (m = 0; m < DIM; m++)
275 if (!opts->nFreeze[gf][m])
277 fam += gmx::square(f[i][m]);
290 for (i = start; i < end; i++)
302 if (la_max >= 0 && DOMAINDECOMP(cr))
304 a_max = cr->dd->globalAtomIndices[la_max];
312 snew(sum, 2*cr->nnodes+1);
313 sum[2*cr->nodeid] = fmax2;
314 sum[2*cr->nodeid+1] = a_max;
315 sum[2*cr->nnodes] = fnorm2;
316 gmx_sumd(2*cr->nnodes+1, sum, cr);
317 fnorm2 = sum[2*cr->nnodes];
318 /* Determine the global maximum */
319 for (i = 0; i < cr->nnodes; i++)
321 if (sum[2*i] > fmax2)
324 a_max = gmx::roundToInt(sum[2*i+1]);
332 *fnorm = sqrt(fnorm2);
344 //! Compute the norm of the force
345 static void get_state_f_norm_max(const t_commrec *cr,
346 t_grpopts *opts, t_mdatoms *mdatoms,
349 get_f_norm_max(cr, opts, mdatoms, ems->f.rvec_array(),
350 &ems->fnorm, &ems->fmax, &ems->a_fmax);
353 //! Initialize the energy minimization
354 static void init_em(FILE *fplog,
355 const gmx::MDLogger &mdlog,
359 gmx::ImdSession *imdSession,
361 t_state *state_global, gmx_mtop_t *top_global,
362 em_state_t *ems, gmx_localtop_t *top,
365 t_graph **graph, gmx::MDAtoms *mdAtoms, gmx_global_stat_t *gstat,
366 gmx_vsite_t *vsite, gmx::Constraints *constr, gmx_shellfc_t **shellfc)
372 fprintf(fplog, "Initiating %s\n", title);
377 state_global->ngtc = 0;
379 initialize_lambdas(fplog, *ir, MASTER(cr), &(state_global->fep_state), state_global->lambda, nullptr);
383 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
385 *shellfc = init_shell_flexcon(stdout,
387 constr ? constr->numFlexibleConstraints() : 0,
393 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
395 /* With energy minimization, shells and flexible constraints are
396 * automatically minimized when treated like normal DOFS.
398 if (shellfc != nullptr)
404 auto mdatoms = mdAtoms->mdatoms();
405 if (DOMAINDECOMP(cr))
407 top->useInDomainDecomp_ = true;
408 dd_init_local_top(*top_global, top);
410 dd_init_local_state(cr->dd, state_global, &ems->s);
412 /* Distribute the charge groups over the nodes from the master node */
413 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
414 state_global, *top_global, ir, imdSession, pull_work,
415 &ems->s, &ems->f, mdAtoms, top,
417 nrnb, nullptr, FALSE);
418 dd_store_state(cr->dd, &ems->s);
424 state_change_natoms(state_global, state_global->natoms);
425 /* Just copy the state */
426 ems->s = *state_global;
427 state_change_natoms(&ems->s, ems->s.natoms);
428 ems->f.resizeWithPadding(ems->s.natoms);
430 mdAlgorithmsSetupAtomData(cr, ir, *top_global, top, fr,
432 constr, vsite, shellfc ? *shellfc : nullptr);
436 set_vsite_top(vsite, top, mdatoms);
440 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
444 // TODO how should this cross-module support dependency be managed?
445 if (ir->eConstrAlg == econtSHAKE &&
446 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
448 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
449 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
452 if (!ir->bContinuation)
454 /* Constrain the starting coordinates */
456 constr->apply(TRUE, TRUE,
458 ems->s.x.rvec_array(),
459 ems->s.x.rvec_array(),
462 ems->s.lambda[efptFEP], &dvdl_constr,
463 nullptr, nullptr, gmx::ConstraintVariable::Positions);
469 *gstat = global_stat_init(ir);
476 calc_shifts(ems->s.box, fr->shift_vec);
479 //! Finalize the minimization
480 static void finish_em(const t_commrec *cr, gmx_mdoutf_t outf,
481 gmx_walltime_accounting_t walltime_accounting,
482 gmx_wallcycle_t wcycle)
484 if (!thisRankHasDuty(cr, DUTY_PME))
486 /* Tell the PME only node to finish */
487 gmx_pme_send_finish(cr);
492 em_time_end(walltime_accounting, wcycle);
495 //! Swap two different EM states during minimization
496 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
505 //! Save the EM trajectory
506 static void write_em_traj(FILE *fplog, const t_commrec *cr,
508 gmx_bool bX, gmx_bool bF, const char *confout,
509 gmx_mtop_t *top_global,
510 t_inputrec *ir, int64_t step,
512 t_state *state_global,
513 ObservablesHistory *observablesHistory)
519 mdof_flags |= MDOF_X;
523 mdof_flags |= MDOF_F;
526 /* If we want IMD output, set appropriate MDOF flag */
529 mdof_flags |= MDOF_IMD;
532 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
533 top_global->natoms, step, static_cast<double>(step),
534 &state->s, state_global, observablesHistory,
537 if (confout != nullptr)
539 if (DOMAINDECOMP(cr))
541 /* If bX=true, x was collected to state_global in the call above */
544 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
545 dd_collect_vec(cr->dd, &state->s, state->s.x, globalXRef);
550 /* Copy the local state pointer */
551 state_global = &state->s;
556 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
558 /* Make molecules whole only for confout writing */
559 do_pbc_mtop(ir->ePBC, state->s.box, top_global,
560 state_global->x.rvec_array());
563 write_sto_conf_mtop(confout,
564 *top_global->name, top_global,
565 state_global->x.rvec_array(), nullptr, ir->ePBC, state->s.box);
570 //! \brief Do one minimization step
572 // \returns true when the step succeeded, false when a constraint error occurred
573 static bool do_em_step(const t_commrec *cr,
574 t_inputrec *ir, t_mdatoms *md,
575 em_state_t *ems1, real a, const PaddedHostVector<gmx::RVec> *force,
577 gmx::Constraints *constr,
584 int nthreads gmx_unused;
586 bool validStep = true;
591 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
593 gmx_incons("state mismatch in do_em_step");
596 s2->flags = s1->flags;
598 if (s2->natoms != s1->natoms)
600 state_change_natoms(s2, s1->natoms);
601 ems2->f.resizeWithPadding(s2->natoms);
603 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
605 s2->cg_gl.resize(s1->cg_gl.size());
608 copy_mat(s1->box, s2->box);
609 /* Copy free energy state */
610 s2->lambda = s1->lambda;
611 copy_mat(s1->box, s2->box);
616 nthreads = gmx_omp_nthreads_get(emntUpdate);
617 #pragma omp parallel num_threads(nthreads)
619 const rvec *x1 = s1->x.rvec_array();
620 rvec *x2 = s2->x.rvec_array();
621 const rvec *f = force->rvec_array();
624 #pragma omp for schedule(static) nowait
625 for (int i = start; i < end; i++)
633 for (int m = 0; m < DIM; m++)
635 if (ir->opts.nFreeze[gf][m])
641 x2[i][m] = x1[i][m] + a*f[i][m];
645 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
648 if (s2->flags & (1<<estCGP))
650 /* Copy the CG p vector */
651 const rvec *p1 = s1->cg_p.rvec_array();
652 rvec *p2 = s2->cg_p.rvec_array();
653 #pragma omp for schedule(static) nowait
654 for (int i = start; i < end; i++)
656 // Trivial OpenMP block that does not throw
657 copy_rvec(p1[i], p2[i]);
661 if (DOMAINDECOMP(cr))
663 s2->ddp_count = s1->ddp_count;
665 /* OpenMP does not supported unsigned loop variables */
666 #pragma omp for schedule(static) nowait
667 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
669 s2->cg_gl[i] = s1->cg_gl[i];
671 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
679 constr->apply(TRUE, TRUE,
681 s1->x.rvec_array(), s2->x.rvec_array(),
683 s2->lambda[efptBONDED], &dvdl_constr,
684 nullptr, nullptr, gmx::ConstraintVariable::Positions);
688 /* This global reduction will affect performance at high
689 * parallelization, but we can not really avoid it.
690 * But usually EM is not run at high parallelization.
692 int reductionBuffer = static_cast<int>(!validStep);
693 gmx_sumi(1, &reductionBuffer, cr);
694 validStep = (reductionBuffer == 0);
697 // We should move this check to the different minimizers
698 if (!validStep && ir->eI != eiSteep)
700 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
701 EI(ir->eI), EI(eiSteep), EI(ir->eI));
708 //! Prepare EM for using domain decomposition parallellization
709 static void em_dd_partition_system(FILE *fplog,
710 const gmx::MDLogger &mdlog,
711 int step, const t_commrec *cr,
712 gmx_mtop_t *top_global, t_inputrec *ir,
713 gmx::ImdSession *imdSession,
715 em_state_t *ems, gmx_localtop_t *top,
716 gmx::MDAtoms *mdAtoms, t_forcerec *fr,
717 gmx_vsite_t *vsite, gmx::Constraints *constr,
718 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
720 /* Repartition the domain decomposition */
721 dd_partition_system(fplog, mdlog, step, cr, FALSE, 1,
722 nullptr, *top_global, ir, imdSession, pull_work,
724 mdAtoms, top, fr, vsite, constr,
725 nrnb, wcycle, FALSE);
726 dd_store_state(cr->dd, &ems->s);
732 /*! \brief Class to handle the work of setting and doing an energy evaluation.
734 * This class is a mere aggregate of parameters to pass to evaluate an
735 * energy, so that future changes to names and types of them consume
736 * less time when refactoring other code.
738 * Aggregate initialization is used, for which the chief risk is that
739 * if a member is added at the end and not all initializer lists are
740 * updated, then the member will be value initialized, which will
741 * typically mean initialization to zero.
743 * Use a braced initializer list to construct one of these. */
744 class EnergyEvaluator
747 /*! \brief Evaluates an energy on the state in \c ems.
749 * \todo In practice, the same objects mu_tot, vir, and pres
750 * are always passed to this function, so we would rather have
751 * them as data members. However, their C-array types are
752 * unsuited for aggregate initialization. When the types
753 * improve, the call signature of this method can be reduced.
755 void run(em_state_t *ems, rvec mu_tot,
756 tensor vir, tensor pres,
757 int64_t count, gmx_bool bFirst);
758 //! Handles logging (deprecated).
761 const gmx::MDLogger &mdlog;
762 //! Handles communication.
764 //! Coordinates multi-simulations.
765 const gmx_multisim_t *ms;
766 //! Holds the simulation topology.
767 gmx_mtop_t *top_global;
768 //! Holds the domain topology.
770 //! User input options.
771 t_inputrec *inputrec;
772 //! The Interactive Molecular Dynamics session.
773 gmx::ImdSession *imdSession;
774 //! The pull work object.
776 //! Manages flop accounting.
778 //! Manages wall cycle accounting.
779 gmx_wallcycle_t wcycle;
780 //! Coordinates global reduction.
781 gmx_global_stat_t gstat;
782 //! Handles virtual sites.
784 //! Handles constraints.
785 gmx::Constraints *constr;
786 //! Handles strange things.
788 //! Molecular graph for SHAKE.
790 //! Per-atom data for this domain.
791 gmx::MDAtoms *mdAtoms;
792 //! Handles how to calculate the forces.
794 //! Schedule of force-calculation work each step for this task.
795 MdrunScheduleWorkload *runScheduleWork;
796 //! Stores the computed energies.
797 gmx_enerdata_t *enerd;
801 EnergyEvaluator::run(em_state_t *ems, rvec mu_tot,
802 tensor vir, tensor pres,
803 int64_t count, gmx_bool bFirst)
807 tensor force_vir, shake_vir, ekin;
811 /* Set the time to the initial time, the time does not change during EM */
812 t = inputrec->init_t;
815 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
817 /* This is the first state or an old state used before the last ns */
823 if (inputrec->nstlist > 0)
831 construct_vsites(vsite, ems->s.x.rvec_array(), 1, nullptr,
832 top->idef.iparams, top->idef.il,
833 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
836 if (DOMAINDECOMP(cr) && bNS)
838 /* Repartition the domain decomposition */
839 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession,
841 ems, top, mdAtoms, fr, vsite, constr,
845 /* Calc force & energy on new trial position */
846 /* do_force always puts the charge groups in the box and shifts again
847 * We do not unshift, so molecules are always whole in congrad.c
849 do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession,
851 count, nrnb, wcycle, top,
852 ems->s.box, ems->s.x.arrayRefWithPadding(), &ems->s.hist,
853 ems->f.arrayRefWithPadding(), force_vir, mdAtoms->mdatoms(), enerd, fcd,
854 ems->s.lambda, graph, fr, runScheduleWork, vsite, mu_tot, t, nullptr,
855 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
856 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
857 (bNS ? GMX_FORCE_NS : 0),
858 DDBalanceRegionHandler(cr));
860 /* Clear the unused shake virial and pressure */
861 clear_mat(shake_vir);
864 /* Communicate stuff when parallel */
865 if (PAR(cr) && inputrec->eI != eiNM)
867 wallcycle_start(wcycle, ewcMoveE);
869 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
870 inputrec, nullptr, nullptr, nullptr, 1, &terminate,
876 wallcycle_stop(wcycle, ewcMoveE);
879 if (fr->dispersionCorrection)
881 /* Calculate long range corrections to pressure and energy */
882 const DispersionCorrection::Correction correction =
883 fr->dispersionCorrection->calculate(ems->s.box, ems->s.lambda[efptVDW]);
885 enerd->term[F_DISPCORR] = correction.energy;
886 enerd->term[F_EPOT] += correction.energy;
887 enerd->term[F_PRES] += correction.pressure;
888 enerd->term[F_DVDL] += correction.dvdl;
892 enerd->term[F_DISPCORR] = 0;
895 ems->epot = enerd->term[F_EPOT];
899 /* Project out the constraint components of the force */
901 rvec *f_rvec = ems->f.rvec_array();
902 constr->apply(FALSE, FALSE,
904 ems->s.x.rvec_array(), f_rvec, f_rvec,
906 ems->s.lambda[efptBONDED], &dvdl_constr,
907 nullptr, &shake_vir, gmx::ConstraintVariable::ForceDispl);
908 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
909 m_add(force_vir, shake_vir, vir);
913 copy_mat(force_vir, vir);
917 enerd->term[F_PRES] =
918 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
920 sum_dhdl(enerd, ems->s.lambda, *inputrec->fepvals);
922 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
924 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
930 //! Parallel utility summing energies and forces
931 static double reorder_partsum(const t_commrec *cr, t_grpopts *opts,
932 gmx_mtop_t *top_global,
933 em_state_t *s_min, em_state_t *s_b)
936 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
941 fprintf(debug, "Doing reorder_partsum\n");
944 const rvec *fm = s_min->f.rvec_array();
945 const rvec *fb = s_b->f.rvec_array();
947 cgs_gl = dd_charge_groups_global(cr->dd);
948 index = cgs_gl->index;
950 /* Collect fm in a global vector fmg.
951 * This conflicts with the spirit of domain decomposition,
952 * but to fully optimize this a much more complicated algorithm is required.
955 snew(fmg, top_global->natoms);
957 ncg = s_min->s.cg_gl.size();
958 cg_gl = s_min->s.cg_gl.data();
960 for (c = 0; c < ncg; c++)
965 for (a = a0; a < a1; a++)
967 copy_rvec(fm[i], fmg[a]);
971 gmx_sum(top_global->natoms*3, fmg[0], cr);
973 /* Now we will determine the part of the sum for the cgs in state s_b */
974 ncg = s_b->s.cg_gl.size();
975 cg_gl = s_b->s.cg_gl.data();
979 gmx::ArrayRef<unsigned char> grpnrFREEZE = top_global->groups.groupNumbers[SimulationAtomGroupType::Freeze];
980 for (c = 0; c < ncg; c++)
985 for (a = a0; a < a1; a++)
987 if (!grpnrFREEZE.empty())
991 for (m = 0; m < DIM; m++)
993 if (!opts->nFreeze[gf][m])
995 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
1007 //! Print some stuff, like beta, whatever that means.
1008 static real pr_beta(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
1009 gmx_mtop_t *top_global,
1010 em_state_t *s_min, em_state_t *s_b)
1014 /* This is just the classical Polak-Ribiere calculation of beta;
1015 * it looks a bit complicated since we take freeze groups into account,
1016 * and might have to sum it in parallel runs.
1019 if (!DOMAINDECOMP(cr) ||
1020 (s_min->s.ddp_count == cr->dd->ddp_count &&
1021 s_b->s.ddp_count == cr->dd->ddp_count))
1023 const rvec *fm = s_min->f.rvec_array();
1024 const rvec *fb = s_b->f.rvec_array();
1027 /* This part of code can be incorrect with DD,
1028 * since the atom ordering in s_b and s_min might differ.
1030 for (int i = 0; i < mdatoms->homenr; i++)
1032 if (mdatoms->cFREEZE)
1034 gf = mdatoms->cFREEZE[i];
1036 for (int m = 0; m < DIM; m++)
1038 if (!opts->nFreeze[gf][m])
1040 sum += (fb[i][m] - fm[i][m])*fb[i][m];
1047 /* We need to reorder cgs while summing */
1048 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1052 gmx_sumd(1, &sum, cr);
1055 return sum/gmx::square(s_min->fnorm);
1062 LegacySimulator::do_cg()
1064 const char *CG = "Polak-Ribiere Conjugate Gradients";
1067 gmx_global_stat_t gstat;
1069 double tmp, minstep;
1071 real a, b, c, beta = 0.0;
1074 gmx_bool converged, foundlower;
1076 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1078 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1079 int m, step, nminstep;
1080 auto mdatoms = mdAtoms->mdatoms();
1082 GMX_LOG(mdlog.info).asParagraph().
1083 appendText("Note that activating conjugate gradient energy minimization via the "
1084 "integrator .mdp option and the command gmx mdrun may "
1085 "be available in a different form in a future version of GROMACS, "
1086 "e.g. gmx minimize and an .mdp option.");
1092 // In CG, the state is extended with a search direction
1093 state_global->flags |= (1<<estCGP);
1095 // Ensure the extra per-atom state array gets allocated
1096 state_change_natoms(state_global, state_global->natoms);
1098 // Initialize the search direction to zero
1099 for (RVec &cg_p : state_global->cg_p)
1105 /* Create 4 states on the stack and extract pointers that we will swap */
1106 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1107 em_state_t *s_min = &s0;
1108 em_state_t *s_a = &s1;
1109 em_state_t *s_b = &s2;
1110 em_state_t *s_c = &s3;
1112 /* Init em and store the local state in s_min */
1113 init_em(fplog, mdlog, CG, cr, inputrec, imdSession,
1115 state_global, top_global, s_min, &top,
1116 nrnb, fr, &graph, mdAtoms, &gstat,
1117 vsite, constr, nullptr);
1118 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1119 StartingBehavior::NewSimulation);
1120 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr, false, mdModulesNotifier);
1122 /* Print to log file */
1123 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1125 /* Max number of steps */
1126 number_steps = inputrec->nsteps;
1130 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1134 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1137 EnergyEvaluator energyEvaluator {
1138 fplog, mdlog, cr, ms,
1140 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1141 vsite, constr, fcd, graph,
1142 mdAtoms, fr, runScheduleWork, enerd
1144 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1145 /* do_force always puts the charge groups in the box and shifts again
1146 * We do not unshift, so molecules are always whole in congrad.c
1148 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1152 /* Copy stuff to the energy bin for easy printing etc. */
1153 matrix nullBox = {};
1154 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1155 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1156 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1158 energyOutput.printHeader(fplog, step, step);
1159 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE,
1164 /* Estimate/guess the initial stepsize */
1165 stepsize = inputrec->em_stepsize/s_min->fnorm;
1169 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1170 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1171 s_min->fmax, s_min->a_fmax+1);
1172 fprintf(stderr, " F-Norm = %12.5e\n",
1173 s_min->fnorm/sqrtNumAtoms);
1174 fprintf(stderr, "\n");
1175 /* and copy to the log file too... */
1176 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1177 s_min->fmax, s_min->a_fmax+1);
1178 fprintf(fplog, " F-Norm = %12.5e\n",
1179 s_min->fnorm/sqrtNumAtoms);
1180 fprintf(fplog, "\n");
1182 /* Start the loop over CG steps.
1183 * Each successful step is counted, and we continue until
1184 * we either converge or reach the max number of steps.
1187 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1190 /* start taking steps in a new direction
1191 * First time we enter the routine, beta=0, and the direction is
1192 * simply the negative gradient.
1195 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1196 rvec *pm = s_min->s.cg_p.rvec_array();
1197 const rvec *sfm = s_min->f.rvec_array();
1200 for (int i = 0; i < mdatoms->homenr; i++)
1202 if (mdatoms->cFREEZE)
1204 gf = mdatoms->cFREEZE[i];
1206 for (m = 0; m < DIM; m++)
1208 if (!inputrec->opts.nFreeze[gf][m])
1210 pm[i][m] = sfm[i][m] + beta*pm[i][m];
1211 gpa -= pm[i][m]*sfm[i][m];
1212 /* f is negative gradient, thus the sign */
1221 /* Sum the gradient along the line across CPUs */
1224 gmx_sumd(1, &gpa, cr);
1227 /* Calculate the norm of the search vector */
1228 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1230 /* Just in case stepsize reaches zero due to numerical precision... */
1233 stepsize = inputrec->em_stepsize/pnorm;
1237 * Double check the value of the derivative in the search direction.
1238 * If it is positive it must be due to the old information in the
1239 * CG formula, so just remove that and start over with beta=0.
1240 * This corresponds to a steepest descent step.
1245 step--; /* Don't count this step since we are restarting */
1246 continue; /* Go back to the beginning of the big for-loop */
1249 /* Calculate minimum allowed stepsize, before the average (norm)
1250 * relative change in coordinate is smaller than precision
1253 auto s_min_x = makeArrayRef(s_min->s.x);
1254 for (int i = 0; i < mdatoms->homenr; i++)
1256 for (m = 0; m < DIM; m++)
1258 tmp = fabs(s_min_x[i][m]);
1267 /* Add up from all CPUs */
1270 gmx_sumd(1, &minstep, cr);
1273 minstep = GMX_REAL_EPS/sqrt(minstep/(3*top_global->natoms));
1275 if (stepsize < minstep)
1281 /* Write coordinates if necessary */
1282 do_x = do_per_step(step, inputrec->nstxout);
1283 do_f = do_per_step(step, inputrec->nstfout);
1285 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1286 top_global, inputrec, step,
1287 s_min, state_global, observablesHistory);
1289 /* Take a step downhill.
1290 * In theory, we should minimize the function along this direction.
1291 * That is quite possible, but it turns out to take 5-10 function evaluations
1292 * for each line. However, we dont really need to find the exact minimum -
1293 * it is much better to start a new CG step in a modified direction as soon
1294 * as we are close to it. This will save a lot of energy evaluations.
1296 * In practice, we just try to take a single step.
1297 * If it worked (i.e. lowered the energy), we increase the stepsize but
1298 * the continue straight to the next CG step without trying to find any minimum.
1299 * If it didn't work (higher energy), there must be a minimum somewhere between
1300 * the old position and the new one.
1302 * Due to the finite numerical accuracy, it turns out that it is a good idea
1303 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1304 * This leads to lower final energies in the tests I've done. / Erik
1306 s_a->epot = s_min->epot;
1308 c = a + stepsize; /* reference position along line is zero */
1310 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1312 em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec, imdSession,
1314 s_min, &top, mdAtoms, fr, vsite, constr,
1318 /* Take a trial step (new coords in s_c) */
1319 do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c,
1323 /* Calculate energy for the trial step */
1324 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1326 /* Calc derivative along line */
1327 const rvec *pc = s_c->s.cg_p.rvec_array();
1328 const rvec *sfc = s_c->f.rvec_array();
1330 for (int i = 0; i < mdatoms->homenr; i++)
1332 for (m = 0; m < DIM; m++)
1334 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1337 /* Sum the gradient along the line across CPUs */
1340 gmx_sumd(1, &gpc, cr);
1343 /* This is the max amount of increase in energy we tolerate */
1344 tmp = std::sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1346 /* Accept the step if the energy is lower, or if it is not significantly higher
1347 * and the line derivative is still negative.
1349 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1352 /* Great, we found a better energy. Increase step for next iteration
1353 * if we are still going down, decrease it otherwise
1357 stepsize *= 1.618034; /* The golden section */
1361 stepsize *= 0.618034; /* 1/golden section */
1366 /* New energy is the same or higher. We will have to do some work
1367 * to find a smaller value in the interval. Take smaller step next time!
1370 stepsize *= 0.618034;
1376 /* OK, if we didn't find a lower value we will have to locate one now - there must
1377 * be one in the interval [a=0,c].
1378 * The same thing is valid here, though: Don't spend dozens of iterations to find
1379 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1380 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1382 * I also have a safeguard for potentially really pathological functions so we never
1383 * take more than 20 steps before we give up ...
1385 * If we already found a lower value we just skip this step and continue to the update.
1394 /* Select a new trial point.
1395 * If the derivatives at points a & c have different sign we interpolate to zero,
1396 * otherwise just do a bisection.
1398 if (gpa < 0 && gpc > 0)
1400 b = a + gpa*(a-c)/(gpc-gpa);
1407 /* safeguard if interpolation close to machine accuracy causes errors:
1408 * never go outside the interval
1410 if (b <= a || b >= c)
1415 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1417 /* Reload the old state */
1418 em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec, imdSession,
1420 s_min, &top, mdAtoms, fr, vsite, constr,
1424 /* Take a trial step to this new point - new coords in s_b */
1425 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b,
1429 /* Calculate energy for the trial step */
1430 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1432 /* p does not change within a step, but since the domain decomposition
1433 * might change, we have to use cg_p of s_b here.
1435 const rvec *pb = s_b->s.cg_p.rvec_array();
1436 const rvec *sfb = s_b->f.rvec_array();
1438 for (int i = 0; i < mdatoms->homenr; i++)
1440 for (m = 0; m < DIM; m++)
1442 gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1445 /* Sum the gradient along the line across CPUs */
1448 gmx_sumd(1, &gpb, cr);
1453 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1454 s_a->epot, s_b->epot, s_c->epot, gpb);
1457 epot_repl = s_b->epot;
1459 /* Keep one of the intervals based on the value of the derivative at the new point */
1462 /* Replace c endpoint with b */
1463 swap_em_state(&s_b, &s_c);
1469 /* Replace a endpoint with b */
1470 swap_em_state(&s_b, &s_a);
1476 * Stop search as soon as we find a value smaller than the endpoints.
1477 * Never run more than 20 steps, no matter what.
1481 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1484 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1487 /* OK. We couldn't find a significantly lower energy.
1488 * If beta==0 this was steepest descent, and then we give up.
1489 * If not, set beta=0 and restart with steepest descent before quitting.
1499 /* Reset memory before giving up */
1505 /* Select min energy state of A & C, put the best in B.
1507 if (s_c->epot < s_a->epot)
1511 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1512 s_c->epot, s_a->epot);
1514 swap_em_state(&s_b, &s_c);
1521 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1522 s_a->epot, s_c->epot);
1524 swap_em_state(&s_b, &s_a);
1533 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1536 swap_em_state(&s_b, &s_c);
1540 /* new search direction */
1541 /* beta = 0 means forget all memory and restart with steepest descents. */
1542 if (nstcg && ((step % nstcg) == 0))
1548 /* s_min->fnorm cannot be zero, because then we would have converged
1552 /* Polak-Ribiere update.
1553 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1555 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1557 /* Limit beta to prevent oscillations */
1558 if (fabs(beta) > 5.0)
1564 /* update positions */
1565 swap_em_state(&s_min, &s_b);
1568 /* Print it if necessary */
1571 if (mdrunOptions.verbose)
1573 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1574 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1575 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1576 s_min->fmax, s_min->a_fmax+1);
1579 /* Store the new (lower) energies */
1580 matrix nullBox = {};
1581 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1582 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1583 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1585 do_log = do_per_step(step, inputrec->nstlog);
1586 do_ene = do_per_step(step, inputrec->nstenergy);
1588 imdSession->fillEnergyRecord(step, TRUE);
1592 energyOutput.printHeader(fplog, step, step);
1594 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1595 do_log ? fplog : nullptr, step, step,
1599 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1600 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0))
1602 imdSession->sendPositionsAndEnergies();
1605 /* Stop when the maximum force lies below tolerance.
1606 * If we have reached machine precision, converged is already set to true.
1608 converged = converged || (s_min->fmax < inputrec->em_tol);
1610 } /* End of the loop */
1614 step--; /* we never took that last step in this case */
1617 if (s_min->fmax > inputrec->em_tol)
1621 warn_step(fplog, inputrec->em_tol, s_min->fmax,
1622 step-1 == number_steps, FALSE);
1629 /* If we printed energy and/or logfile last step (which was the last step)
1630 * we don't have to do it again, but otherwise print the final values.
1634 /* Write final value to log since we didn't do anything the last step */
1635 energyOutput.printHeader(fplog, step, step);
1637 if (!do_ene || !do_log)
1639 /* Write final energy file entries */
1640 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1641 !do_log ? fplog : nullptr, step, step,
1646 /* Print some stuff... */
1649 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1653 * For accurate normal mode calculation it is imperative that we
1654 * store the last conformation into the full precision binary trajectory.
1656 * However, we should only do it if we did NOT already write this step
1657 * above (which we did if do_x or do_f was true).
1659 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1660 * in the trajectory with the same step number.
1662 do_x = !do_per_step(step, inputrec->nstxout);
1663 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1665 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1666 top_global, inputrec, step,
1667 s_min, state_global, observablesHistory);
1672 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1673 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1674 s_min, sqrtNumAtoms);
1675 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1676 s_min, sqrtNumAtoms);
1678 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1681 finish_em(cr, outf, walltime_accounting, wcycle);
1683 /* To print the actual number of steps we needed somewhere */
1684 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1689 LegacySimulator::do_lbfgs()
1691 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1694 gmx_global_stat_t gstat;
1696 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1697 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1698 real *rho, *alpha, *p, *s, **dx, **dg;
1699 real a, b, c, maxdelta, delta;
1701 real dgdx, dgdg, sq, yr, beta;
1704 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1706 int start, end, number_steps;
1707 int i, k, m, n, gf, step;
1709 auto mdatoms = mdAtoms->mdatoms();
1711 GMX_LOG(mdlog.info).asParagraph().
1712 appendText("Note that activating L-BFGS energy minimization via the "
1713 "integrator .mdp option and the command gmx mdrun may "
1714 "be available in a different form in a future version of GROMACS, "
1715 "e.g. gmx minimize and an .mdp option.");
1719 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1722 if (nullptr != constr)
1724 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).");
1727 n = 3*state_global->natoms;
1728 nmaxcorr = inputrec->nbfgscorr;
1733 snew(rho, nmaxcorr);
1734 snew(alpha, nmaxcorr);
1737 for (i = 0; i < nmaxcorr; i++)
1743 for (i = 0; i < nmaxcorr; i++)
1752 init_em(fplog, mdlog, LBFGS, cr, inputrec, imdSession,
1754 state_global, top_global, &ems, &top,
1755 nrnb, fr, &graph, mdAtoms, &gstat,
1756 vsite, constr, nullptr);
1757 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1758 StartingBehavior::NewSimulation);
1759 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr, false, mdModulesNotifier);
1762 end = mdatoms->homenr;
1764 /* We need 4 working states */
1765 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1766 em_state_t *sa = &s0;
1767 em_state_t *sb = &s1;
1768 em_state_t *sc = &s2;
1769 em_state_t *last = &s3;
1770 /* Initialize by copying the state from ems (we could skip x and f here) */
1775 /* Print to log file */
1776 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1778 do_log = do_ene = do_x = do_f = TRUE;
1780 /* Max number of steps */
1781 number_steps = inputrec->nsteps;
1783 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1785 for (i = start; i < end; i++)
1787 if (mdatoms->cFREEZE)
1789 gf = mdatoms->cFREEZE[i];
1791 for (m = 0; m < DIM; m++)
1793 frozen[3*i+m] = (inputrec->opts.nFreeze[gf][m] != 0);
1798 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1802 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1807 construct_vsites(vsite, state_global->x.rvec_array(), 1, nullptr,
1808 top.idef.iparams, top.idef.il,
1809 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1812 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1813 /* do_force always puts the charge groups in the box and shifts again
1814 * We do not unshift, so molecules are always whole
1817 EnergyEvaluator energyEvaluator {
1818 fplog, mdlog, cr, ms,
1820 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1821 vsite, constr, fcd, graph,
1822 mdAtoms, fr, runScheduleWork, enerd
1824 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1828 /* Copy stuff to the energy bin for easy printing etc. */
1829 matrix nullBox = {};
1830 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1831 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1832 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1834 energyOutput.printHeader(fplog, step, step);
1835 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE,
1840 /* Set the initial step.
1841 * since it will be multiplied by the non-normalized search direction
1842 * vector (force vector the first time), we scale it by the
1843 * norm of the force.
1848 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1849 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1850 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1851 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1852 fprintf(stderr, "\n");
1853 /* and copy to the log file too... */
1854 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1855 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1856 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1857 fprintf(fplog, "\n");
1860 // Point is an index to the memory of search directions, where 0 is the first one.
1863 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1864 real *fInit = static_cast<real *>(ems.f.rvec_array()[0]);
1865 for (i = 0; i < n; i++)
1869 dx[point][i] = fInit[i]; /* Initial search direction */
1877 // Stepsize will be modified during the search, and actually it is not critical
1878 // (the main efficiency in the algorithm comes from changing directions), but
1879 // we still need an initial value, so estimate it as the inverse of the norm
1880 // so we take small steps where the potential fluctuates a lot.
1881 stepsize = 1.0/ems.fnorm;
1883 /* Start the loop over BFGS steps.
1884 * Each successful step is counted, and we continue until
1885 * we either converge or reach the max number of steps.
1890 /* Set the gradient from the force */
1892 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1895 /* Write coordinates if necessary */
1896 do_x = do_per_step(step, inputrec->nstxout);
1897 do_f = do_per_step(step, inputrec->nstfout);
1902 mdof_flags |= MDOF_X;
1907 mdof_flags |= MDOF_F;
1912 mdof_flags |= MDOF_IMD;
1915 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1916 top_global->natoms, step, static_cast<real>(step), &ems.s,
1917 state_global, observablesHistory, ems.f);
1919 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1921 /* make s a pointer to current search direction - point=0 first time we get here */
1924 real *xx = static_cast<real *>(ems.s.x.rvec_array()[0]);
1925 real *ff = static_cast<real *>(ems.f.rvec_array()[0]);
1927 // calculate line gradient in position A
1928 for (gpa = 0, i = 0; i < n; i++)
1933 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1934 * relative change in coordinate is smaller than precision
1936 for (minstep = 0, i = 0; i < n; i++)
1946 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1948 if (stepsize < minstep)
1954 // Before taking any steps along the line, store the old position
1956 real *lastx = static_cast<real *>(last->s.x.data()[0]);
1957 real *lastf = static_cast<real *>(last->f.data()[0]);
1962 /* Take a step downhill.
1963 * In theory, we should find the actual minimum of the function in this
1964 * direction, somewhere along the line.
1965 * That is quite possible, but it turns out to take 5-10 function evaluations
1966 * for each line. However, we dont really need to find the exact minimum -
1967 * it is much better to start a new BFGS step in a modified direction as soon
1968 * as we are close to it. This will save a lot of energy evaluations.
1970 * In practice, we just try to take a single step.
1971 * If it worked (i.e. lowered the energy), we increase the stepsize but
1972 * continue straight to the next BFGS step without trying to find any minimum,
1973 * i.e. we change the search direction too. If the line was smooth, it is
1974 * likely we are in a smooth region, and then it makes sense to take longer
1975 * steps in the modified search direction too.
1977 * If it didn't work (higher energy), there must be a minimum somewhere between
1978 * the old position and the new one. Then we need to start by finding a lower
1979 * value before we change search direction. Since the energy was apparently
1980 * quite rough, we need to decrease the step size.
1982 * Due to the finite numerical accuracy, it turns out that it is a good idea
1983 * to accept a SMALL increase in energy, if the derivative is still downhill.
1984 * This leads to lower final energies in the tests I've done. / Erik
1987 // State "A" is the first position along the line.
1988 // reference position along line is initially zero
1991 // Check stepsize first. We do not allow displacements
1992 // larger than emstep.
1996 // Pick a new position C by adding stepsize to A.
1999 // Calculate what the largest change in any individual coordinate
2000 // would be (translation along line * gradient along line)
2002 for (i = 0; i < n; i++)
2005 if (delta > maxdelta)
2010 // If any displacement is larger than the stepsize limit, reduce the step
2011 if (maxdelta > inputrec->em_stepsize)
2016 while (maxdelta > inputrec->em_stepsize);
2018 // Take a trial step and move the coordinate array xc[] to position C
2019 real *xc = static_cast<real *>(sc->s.x.rvec_array()[0]);
2020 for (i = 0; i < n; i++)
2022 xc[i] = lastx[i] + c*s[i];
2026 // Calculate energy for the trial step in position C
2027 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2029 // Calc line gradient in position C
2030 real *fc = static_cast<real *>(sc->f.rvec_array()[0]);
2031 for (gpc = 0, i = 0; i < n; i++)
2033 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2035 /* Sum the gradient along the line across CPUs */
2038 gmx_sumd(1, &gpc, cr);
2041 // This is the max amount of increase in energy we tolerate.
2042 // By allowing VERY small changes (close to numerical precision) we
2043 // frequently find even better (lower) final energies.
2044 tmp = std::sqrt(GMX_REAL_EPS)*fabs(sa->epot);
2046 // Accept the step if the energy is lower in the new position C (compared to A),
2047 // or if it is not significantly higher and the line derivative is still negative.
2048 foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2049 // If true, great, we found a better energy. We no longer try to alter the
2050 // stepsize, but simply accept this new better position. The we select a new
2051 // search direction instead, which will be much more efficient than continuing
2052 // to take smaller steps along a line. Set fnorm based on the new C position,
2053 // which will be used to update the stepsize to 1/fnorm further down.
2055 // If false, the energy is NOT lower in point C, i.e. it will be the same
2056 // or higher than in point A. In this case it is pointless to move to point C,
2057 // so we will have to do more iterations along the same line to find a smaller
2058 // value in the interval [A=0.0,C].
2059 // Here, A is still 0.0, but that will change when we do a search in the interval
2060 // [0.0,C] below. That search we will do by interpolation or bisection rather
2061 // than with the stepsize, so no need to modify it. For the next search direction
2062 // it will be reset to 1/fnorm anyway.
2066 // OK, if we didn't find a lower value we will have to locate one now - there must
2067 // be one in the interval [a,c].
2068 // The same thing is valid here, though: Don't spend dozens of iterations to find
2069 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2070 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2071 // I also have a safeguard for potentially really pathological functions so we never
2072 // take more than 20 steps before we give up.
2073 // If we already found a lower value we just skip this step and continue to the update.
2078 // Select a new trial point B in the interval [A,C].
2079 // If the derivatives at points a & c have different sign we interpolate to zero,
2080 // otherwise just do a bisection since there might be multiple minima/maxima
2081 // inside the interval.
2082 if (gpa < 0 && gpc > 0)
2084 b = a + gpa*(a-c)/(gpc-gpa);
2091 /* safeguard if interpolation close to machine accuracy causes errors:
2092 * never go outside the interval
2094 if (b <= a || b >= c)
2099 // Take a trial step to point B
2100 real *xb = static_cast<real *>(sb->s.x.rvec_array()[0]);
2101 for (i = 0; i < n; i++)
2103 xb[i] = lastx[i] + b*s[i];
2107 // Calculate energy for the trial step in point B
2108 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2111 // Calculate gradient in point B
2112 real *fb = static_cast<real *>(sb->f.rvec_array()[0]);
2113 for (gpb = 0, i = 0; i < n; i++)
2115 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2118 /* Sum the gradient along the line across CPUs */
2121 gmx_sumd(1, &gpb, cr);
2124 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2125 // at the new point B, and rename the endpoints of this new interval A and C.
2128 /* Replace c endpoint with b */
2130 /* copy state b to c */
2135 /* Replace a endpoint with b */
2137 /* copy state b to a */
2142 * Stop search as soon as we find a value smaller than the endpoints,
2143 * or if the tolerance is below machine precision.
2144 * Never run more than 20 steps, no matter what.
2148 while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2150 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2152 /* OK. We couldn't find a significantly lower energy.
2153 * If ncorr==0 this was steepest descent, and then we give up.
2154 * If not, reset memory to restart as steepest descent before quitting.
2166 /* Search in gradient direction */
2167 for (i = 0; i < n; i++)
2169 dx[point][i] = ff[i];
2171 /* Reset stepsize */
2172 stepsize = 1.0/fnorm;
2177 /* Select min energy state of A & C, put the best in xx/ff/Epot
2179 if (sc->epot < sa->epot)
2201 /* Update the memory information, and calculate a new
2202 * approximation of the inverse hessian
2205 /* Have new data in Epot, xx, ff */
2206 if (ncorr < nmaxcorr)
2211 for (i = 0; i < n; i++)
2213 dg[point][i] = lastf[i]-ff[i];
2214 dx[point][i] *= step_taken;
2219 for (i = 0; i < n; i++)
2221 dgdg += dg[point][i]*dg[point][i];
2222 dgdx += dg[point][i]*dx[point][i];
2227 rho[point] = 1.0/dgdx;
2230 if (point >= nmaxcorr)
2236 for (i = 0; i < n; i++)
2243 /* Recursive update. First go back over the memory points */
2244 for (k = 0; k < ncorr; k++)
2253 for (i = 0; i < n; i++)
2255 sq += dx[cp][i]*p[i];
2258 alpha[cp] = rho[cp]*sq;
2260 for (i = 0; i < n; i++)
2262 p[i] -= alpha[cp]*dg[cp][i];
2266 for (i = 0; i < n; i++)
2271 /* And then go forward again */
2272 for (k = 0; k < ncorr; k++)
2275 for (i = 0; i < n; i++)
2277 yr += p[i]*dg[cp][i];
2281 beta = alpha[cp]-beta;
2283 for (i = 0; i < n; i++)
2285 p[i] += beta*dx[cp][i];
2295 for (i = 0; i < n; i++)
2299 dx[point][i] = p[i];
2307 /* Print it if necessary */
2310 if (mdrunOptions.verbose)
2312 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2313 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2314 step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2317 /* Store the new (lower) energies */
2318 matrix nullBox = {};
2319 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
2320 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2321 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2323 do_log = do_per_step(step, inputrec->nstlog);
2324 do_ene = do_per_step(step, inputrec->nstenergy);
2326 imdSession->fillEnergyRecord(step, TRUE);
2330 energyOutput.printHeader(fplog, step, step);
2332 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2333 do_log ? fplog : nullptr, step, step,
2337 /* Send x and E to IMD client, if bIMD is TRUE. */
2338 if (imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0) && MASTER(cr))
2340 imdSession->sendPositionsAndEnergies();
2343 // Reset stepsize in we are doing more iterations
2346 /* Stop when the maximum force lies below tolerance.
2347 * If we have reached machine precision, converged is already set to true.
2349 converged = converged || (ems.fmax < inputrec->em_tol);
2351 } /* End of the loop */
2355 step--; /* we never took that last step in this case */
2358 if (ems.fmax > inputrec->em_tol)
2362 warn_step(fplog, inputrec->em_tol, ems.fmax,
2363 step-1 == number_steps, FALSE);
2368 /* If we printed energy and/or logfile last step (which was the last step)
2369 * we don't have to do it again, but otherwise print the final values.
2371 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2373 energyOutput.printHeader(fplog, step, step);
2375 if (!do_ene || !do_log) /* Write final energy file entries */
2377 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2378 !do_log ? fplog : nullptr, step, step,
2382 /* Print some stuff... */
2385 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2389 * For accurate normal mode calculation it is imperative that we
2390 * store the last conformation into the full precision binary trajectory.
2392 * However, we should only do it if we did NOT already write this step
2393 * above (which we did if do_x or do_f was true).
2395 do_x = !do_per_step(step, inputrec->nstxout);
2396 do_f = !do_per_step(step, inputrec->nstfout);
2397 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2398 top_global, inputrec, step,
2399 &ems, state_global, observablesHistory);
2403 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2404 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2405 number_steps, &ems, sqrtNumAtoms);
2406 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2407 number_steps, &ems, sqrtNumAtoms);
2409 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2412 finish_em(cr, outf, walltime_accounting, wcycle);
2414 /* To print the actual number of steps we needed somewhere */
2415 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2419 LegacySimulator::do_steep()
2421 const char *SD = "Steepest Descents";
2423 gmx_global_stat_t gstat;
2427 gmx_bool bDone, bAbort, do_x, do_f;
2432 int steps_accepted = 0;
2433 auto mdatoms = mdAtoms->mdatoms();
2435 GMX_LOG(mdlog.info).asParagraph().
2436 appendText("Note that activating steepest-descent energy minimization via the "
2437 "integrator .mdp option and the command gmx mdrun may "
2438 "be available in a different form in a future version of GROMACS, "
2439 "e.g. gmx minimize and an .mdp option.");
2441 /* Create 2 states on the stack and extract pointers that we will swap */
2442 em_state_t s0 {}, s1 {};
2443 em_state_t *s_min = &s0;
2444 em_state_t *s_try = &s1;
2446 /* Init em and store the local state in s_try */
2447 init_em(fplog, mdlog, SD, cr, inputrec, imdSession,
2449 state_global, top_global, s_try, &top,
2450 nrnb, fr, &graph, mdAtoms, &gstat,
2451 vsite, constr, nullptr);
2452 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2453 StartingBehavior::NewSimulation);
2454 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr, false, mdModulesNotifier);
2456 /* Print to log file */
2457 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2459 /* Set variables for stepsize (in nm). This is the largest
2460 * step that we are going to make in any direction.
2462 ustep = inputrec->em_stepsize;
2465 /* Max number of steps */
2466 nsteps = inputrec->nsteps;
2470 /* Print to the screen */
2471 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2475 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2477 EnergyEvaluator energyEvaluator {
2478 fplog, mdlog, cr, ms,
2480 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2481 vsite, constr, fcd, graph,
2482 mdAtoms, fr, runScheduleWork, enerd
2485 /**** HERE STARTS THE LOOP ****
2486 * count is the counter for the number of steps
2487 * bDone will be TRUE when the minimization has converged
2488 * bAbort will be TRUE when nsteps steps have been performed or when
2489 * the stepsize becomes smaller than is reasonable for machine precision
2494 while (!bDone && !bAbort)
2496 bAbort = (nsteps >= 0) && (count == nsteps);
2498 /* set new coordinates, except for first step */
2499 bool validStep = true;
2503 do_em_step(cr, inputrec, mdatoms,
2504 s_min, stepsize, &s_min->f, s_try,
2510 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2514 // Signal constraint error during stepping with energy=inf
2515 s_try->epot = std::numeric_limits<real>::infinity();
2520 energyOutput.printHeader(fplog, count, count);
2525 s_min->epot = s_try->epot;
2528 /* Print it if necessary */
2531 if (mdrunOptions.verbose)
2533 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2534 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2535 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2539 if ( (count == 0) || (s_try->epot < s_min->epot) )
2541 /* Store the new (lower) energies */
2542 matrix nullBox = {};
2543 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(count),
2544 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2545 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2547 imdSession->fillEnergyRecord(count, TRUE);
2549 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2550 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2551 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE,
2553 fplog, count, count,
2559 /* Now if the new energy is smaller than the previous...
2560 * or if this is the first step!
2561 * or if we did random steps!
2564 if ( (count == 0) || (s_try->epot < s_min->epot) )
2568 /* Test whether the convergence criterion is met... */
2569 bDone = (s_try->fmax < inputrec->em_tol);
2571 /* Copy the arrays for force, positions and energy */
2572 /* The 'Min' array always holds the coords and forces of the minimal
2574 swap_em_state(&s_min, &s_try);
2580 /* Write to trn, if necessary */
2581 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2582 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2583 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2584 top_global, inputrec, count,
2585 s_min, state_global, observablesHistory);
2589 /* If energy is not smaller make the step smaller... */
2592 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2594 /* Reload the old state */
2595 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession,
2597 s_min, &top, mdAtoms, fr, vsite, constr,
2602 // If the force is very small after finishing minimization,
2603 // we risk dividing by zero when calculating the step size.
2604 // So we check first if the minimization has stopped before
2605 // trying to obtain a new step size.
2608 /* Determine new step */
2609 stepsize = ustep/s_min->fmax;
2612 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2614 if (count == nsteps || ustep < 1e-12)
2616 if (count == nsteps || ustep < 1e-6)
2621 warn_step(fplog, inputrec->em_tol, s_min->fmax,
2622 count == nsteps, constr != nullptr);
2627 /* Send IMD energies and positions, if bIMD is TRUE. */
2628 if (imdSession->run(count, TRUE, state_global->box,
2629 MASTER(cr) ? state_global->x.rvec_array() : nullptr,
2633 imdSession->sendPositionsAndEnergies();
2637 } /* End of the loop */
2639 /* Print some data... */
2642 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2644 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2645 top_global, inputrec, count,
2646 s_min, state_global, observablesHistory);
2650 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2652 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2653 s_min, sqrtNumAtoms);
2654 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2655 s_min, sqrtNumAtoms);
2658 finish_em(cr, outf, walltime_accounting, wcycle);
2660 /* To print the actual number of steps we needed somewhere */
2661 inputrec->nsteps = count;
2663 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2667 LegacySimulator::do_nm()
2669 const char *NM = "Normal Mode Analysis";
2672 gmx_global_stat_t gstat;
2677 gmx_bool bSparse; /* use sparse matrix storage format */
2679 gmx_sparsematrix_t * sparse_matrix = nullptr;
2680 real * full_matrix = nullptr;
2682 /* added with respect to mdrun */
2684 real der_range = 10.0*std::sqrt(GMX_REAL_EPS);
2686 bool bIsMaster = MASTER(cr);
2687 auto mdatoms = mdAtoms->mdatoms();
2689 GMX_LOG(mdlog.info).asParagraph().
2690 appendText("Note that activating normal-mode analysis via the integrator "
2691 ".mdp option and the command gmx mdrun may "
2692 "be available in a different form in a future version of GROMACS, "
2693 "e.g. gmx normal-modes.");
2695 if (constr != nullptr)
2697 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2700 gmx_shellfc_t *shellfc;
2702 em_state_t state_work {};
2704 /* Init em and store the local state in state_minimum */
2705 init_em(fplog, mdlog, NM, cr, inputrec, imdSession,
2707 state_global, top_global, &state_work, &top,
2708 nrnb, fr, &graph, mdAtoms, &gstat,
2709 vsite, constr, &shellfc);
2710 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2711 StartingBehavior::NewSimulation);
2713 std::vector<int> atom_index = get_atom_index(top_global);
2714 std::vector<gmx::RVec> fneg(atom_index.size(), {0, 0, 0});
2715 snew(dfdx, atom_index.size());
2721 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2722 " which MIGHT not be accurate enough for normal mode analysis.\n"
2723 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2724 " are fairly modest even if you recompile in double precision.\n\n");
2728 /* Check if we can/should use sparse storage format.
2730 * Sparse format is only useful when the Hessian itself is sparse, which it
2731 * will be when we use a cutoff.
2732 * For small systems (n<1000) it is easier to always use full matrix format, though.
2734 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2736 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2739 else if (atom_index.size() < 1000)
2741 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2747 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2751 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2752 sz = DIM*atom_index.size();
2754 fprintf(stderr, "Allocating Hessian memory...\n\n");
2758 sparse_matrix = gmx_sparsematrix_init(sz);
2759 sparse_matrix->compressed_symmetric = TRUE;
2763 snew(full_matrix, sz*sz);
2766 /* Write start time and temperature */
2767 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2769 /* fudge nr of steps to nr of atoms */
2770 inputrec->nsteps = atom_index.size()*2;
2774 fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2775 *(top_global->name), inputrec->nsteps);
2778 nnodes = cr->nnodes;
2780 /* Make evaluate_energy do a single node force calculation */
2782 EnergyEvaluator energyEvaluator {
2783 fplog, mdlog, cr, ms,
2785 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2786 vsite, constr, fcd, graph,
2787 mdAtoms, fr, runScheduleWork, enerd
2789 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2790 cr->nnodes = nnodes;
2792 /* if forces are not small, warn user */
2793 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2795 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2796 if (state_work.fmax > 1.0e-3)
2798 GMX_LOG(mdlog.warning).appendText(
2799 "The force is probably not small enough to "
2800 "ensure that you are at a minimum.\n"
2801 "Be aware that negative eigenvalues may occur\n"
2802 "when the resulting matrix is diagonalized.");
2805 /***********************************************************
2807 * Loop over all pairs in matrix
2809 * do_force called twice. Once with positive and
2810 * once with negative displacement
2812 ************************************************************/
2814 /* Steps are divided one by one over the nodes */
2816 auto state_work_x = makeArrayRef(state_work.s.x);
2817 auto state_work_f = makeArrayRef(state_work.f);
2818 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
2820 size_t atom = atom_index[aid];
2821 for (size_t d = 0; d < DIM; d++)
2824 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2827 x_min = state_work_x[atom][d];
2829 for (unsigned int dx = 0; (dx < 2); dx++)
2833 state_work_x[atom][d] = x_min - der_range;
2837 state_work_x[atom][d] = x_min + der_range;
2840 /* Make evaluate_energy do a single node force calculation */
2844 /* Now is the time to relax the shells */
2845 relax_shell_flexcon(fplog,
2848 mdrunOptions.verbose,
2860 state_work.s.natoms,
2861 state_work.s.x.arrayRefWithPadding(),
2862 state_work.s.v.arrayRefWithPadding(),
2864 state_work.s.lambda,
2866 state_work.f.arrayRefWithPadding(),
2878 DDBalanceRegionHandler(nullptr));
2884 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid*2+dx, FALSE);
2887 cr->nnodes = nnodes;
2891 std::copy(state_work_f.begin(), state_work_f.begin()+atom_index.size(), fneg.begin());
2895 /* x is restored to original */
2896 state_work_x[atom][d] = x_min;
2898 for (size_t j = 0; j < atom_index.size(); j++)
2900 for (size_t k = 0; (k < DIM); k++)
2903 -(state_work_f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2910 #define mpi_type GMX_MPI_REAL
2911 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2912 cr->nodeid, cr->mpi_comm_mygroup);
2917 for (index node = 0; (node < nnodes && aid+node < ssize(atom_index)); node++)
2923 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2924 cr->mpi_comm_mygroup, &stat);
2929 row = (aid + node)*DIM + d;
2931 for (size_t j = 0; j < atom_index.size(); j++)
2933 for (size_t k = 0; k < DIM; k++)
2939 if (col >= row && dfdx[j][k] != 0.0)
2941 gmx_sparsematrix_increment_value(sparse_matrix,
2942 row, col, dfdx[j][k]);
2947 full_matrix[row*sz+col] = dfdx[j][k];
2954 if (mdrunOptions.verbose && fplog)
2959 /* write progress */
2960 if (bIsMaster && mdrunOptions.verbose)
2962 fprintf(stderr, "\rFinished step %d out of %td",
2963 std::min<int>(atom+nnodes, atom_index.size()),
2971 fprintf(stderr, "\n\nWriting Hessian...\n");
2972 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2975 finish_em(cr, outf, walltime_accounting, wcycle);
2977 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);