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 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 * \brief This file defines integrators for energy minimization
42 * \author Berk Hess <hess@kth.se>
43 * \author Erik Lindahl <erik@kth.se>
44 * \ingroup module_mdrun
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/collect.h"
59 #include "gromacs/domdec/dlbtiming.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/mdsetup.h"
63 #include "gromacs/domdec/partition.h"
64 #include "gromacs/ewald/pme_pp.h"
65 #include "gromacs/fileio/confio.h"
66 #include "gromacs/fileio/mtxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/imd/imd.h"
70 #include "gromacs/linearalgebra/sparsematrix.h"
71 #include "gromacs/listed_forces/listed_forces.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/mdlib/constr.h"
75 #include "gromacs/mdlib/coupling.h"
76 #include "gromacs/mdlib/dispersioncorrection.h"
77 #include "gromacs/mdlib/ebin.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/energyoutput.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/gmx_omp_nthreads.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/stat.h"
87 #include "gromacs/mdlib/tgroup.h"
88 #include "gromacs/mdlib/trajectory_writing.h"
89 #include "gromacs/mdlib/update.h"
90 #include "gromacs/mdlib/vsite.h"
91 #include "gromacs/mdrunutility/handlerestart.h"
92 #include "gromacs/mdrunutility/printtime.h"
93 #include "gromacs/mdtypes/commrec.h"
94 #include "gromacs/mdtypes/forcerec.h"
95 #include "gromacs/mdtypes/inputrec.h"
96 #include "gromacs/mdtypes/interaction_const.h"
97 #include "gromacs/mdtypes/md_enums.h"
98 #include "gromacs/mdtypes/mdatom.h"
99 #include "gromacs/mdtypes/mdrunoptions.h"
100 #include "gromacs/mdtypes/state.h"
101 #include "gromacs/pbcutil/pbc.h"
102 #include "gromacs/timing/wallcycle.h"
103 #include "gromacs/timing/walltime_accounting.h"
104 #include "gromacs/topology/mtop_util.h"
105 #include "gromacs/topology/topology.h"
106 #include "gromacs/utility/cstringutil.h"
107 #include "gromacs/utility/exceptions.h"
108 #include "gromacs/utility/fatalerror.h"
109 #include "gromacs/utility/logger.h"
110 #include "gromacs/utility/smalloc.h"
112 #include "legacysimulator.h"
116 using gmx::MdrunScheduleWorkload;
118 using gmx::VirtualSitesHandler;
120 //! Utility structure for manipulating states during EM
121 typedef struct em_state
123 //! Copy of the global state
126 PaddedHostVector<gmx::RVec> f;
129 //! Norm of the force
137 //! Print the EM starting conditions
138 static void print_em_start(FILE* fplog,
140 gmx_walltime_accounting_t walltime_accounting,
141 gmx_wallcycle_t wcycle,
144 walltime_accounting_start_time(walltime_accounting);
145 wallcycle_start(wcycle, ewcRUN);
146 print_start(fplog, cr, walltime_accounting, name);
149 //! Stop counting time for EM
150 static void em_time_end(gmx_walltime_accounting_t walltime_accounting, gmx_wallcycle_t wcycle)
152 wallcycle_stop(wcycle, ewcRUN);
154 walltime_accounting_end_time(walltime_accounting);
157 //! Printing a log file and console header
158 static void sp_header(FILE* out, const char* minimizer, real ftol, int nsteps)
161 fprintf(out, "%s:\n", minimizer);
162 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
163 fprintf(out, " Number of steps = %12d\n", nsteps);
166 //! Print warning message
167 static void warn_step(FILE* fp, real ftol, real fmax, gmx_bool bLastStep, gmx_bool bConstrain)
169 constexpr bool realIsDouble = GMX_DOUBLE;
172 if (!std::isfinite(fmax))
175 "\nEnergy minimization has stopped because the force "
176 "on at least one atom is not finite. This usually means "
177 "atoms are overlapping. Modify the input coordinates to "
178 "remove atom overlap or use soft-core potentials with "
179 "the free energy code to avoid infinite forces.\n%s",
180 !realIsDouble ? "You could also be lucky that switching to double precision "
181 "is sufficient to obtain finite forces.\n"
187 "\nEnergy minimization reached the maximum number "
188 "of steps before the forces reached the requested "
189 "precision Fmax < %g.\n",
195 "\nEnergy minimization has stopped, but the forces have "
196 "not converged to the requested precision Fmax < %g (which "
197 "may not be possible for your system). It stopped "
198 "because the algorithm tried to make a new step whose size "
199 "was too small, or there was no change in the energy since "
200 "last step. Either way, we regard the minimization as "
201 "converged to within the available machine precision, "
202 "given your starting configuration and EM parameters.\n%s%s",
204 !realIsDouble ? "\nDouble precision normally gives you higher accuracy, but "
205 "this is often not needed for preparing to run molecular "
208 bConstrain ? "You might need to increase your constraint accuracy, or turn\n"
209 "off constraints altogether (set constraints = none in mdp file)\n"
213 fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
214 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
217 //! Print message about convergence of the EM
218 static void print_converged(FILE* fp,
224 const em_state_t* ems,
227 char buf[STEPSTRSIZE];
231 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n", alg, ftol, gmx_step_str(count, buf));
233 else if (count < nsteps)
236 "\n%s converged to machine precision in %s steps,\n"
237 "but did not reach the requested Fmax < %g.\n",
238 alg, gmx_step_str(count, buf), ftol);
242 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n", alg, ftol,
243 gmx_step_str(count, buf));
247 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
248 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
249 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm / sqrtNumAtoms);
251 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
252 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
253 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm / sqrtNumAtoms);
257 //! Compute the norm and max of the force array in parallel
258 static void get_f_norm_max(const t_commrec* cr,
268 int la_max, a_max, start, end, i, m, gf;
270 /* This routine finds the largest force and returns it.
271 * On parallel machines the global max is taken.
277 end = mdatoms->homenr;
278 if (mdatoms->cFREEZE)
280 for (i = start; i < end; i++)
282 gf = mdatoms->cFREEZE[i];
284 for (m = 0; m < DIM; m++)
286 if (!opts->nFreeze[gf][m])
288 fam += gmx::square(f[i][m]);
301 for (i = start; i < end; i++)
313 if (la_max >= 0 && DOMAINDECOMP(cr))
315 a_max = cr->dd->globalAtomIndices[la_max];
323 snew(sum, 2 * cr->nnodes + 1);
324 sum[2 * cr->nodeid] = fmax2;
325 sum[2 * cr->nodeid + 1] = a_max;
326 sum[2 * cr->nnodes] = fnorm2;
327 gmx_sumd(2 * cr->nnodes + 1, sum, cr);
328 fnorm2 = sum[2 * cr->nnodes];
329 /* Determine the global maximum */
330 for (i = 0; i < cr->nnodes; i++)
332 if (sum[2 * i] > fmax2)
335 a_max = gmx::roundToInt(sum[2 * i + 1]);
343 *fnorm = sqrt(fnorm2);
355 //! Compute the norm of the force
356 static void get_state_f_norm_max(const t_commrec* cr, t_grpopts* opts, t_mdatoms* mdatoms, em_state_t* ems)
358 get_f_norm_max(cr, opts, mdatoms, ems->f.rvec_array(), &ems->fnorm, &ems->fmax, &ems->a_fmax);
361 //! Initialize the energy minimization
362 static void init_em(FILE* fplog,
363 const gmx::MDLogger& mdlog,
367 gmx::ImdSession* imdSession,
369 t_state* state_global,
370 const gmx_mtop_t* top_global,
375 gmx::MDAtoms* mdAtoms,
376 gmx_global_stat_t* gstat,
377 VirtualSitesHandler* vsite,
378 gmx::Constraints* constr,
379 gmx_shellfc_t** shellfc)
385 fprintf(fplog, "Initiating %s\n", title);
390 state_global->ngtc = 0;
392 initialize_lambdas(fplog, *ir, MASTER(cr), &(state_global->fep_state), state_global->lambda, nullptr);
396 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
398 *shellfc = init_shell_flexcon(stdout, top_global, constr ? constr->numFlexibleConstraints() : 0,
399 ir->nstcalcenergy, DOMAINDECOMP(cr));
403 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI),
404 "This else currently only handles energy minimizers, consider if your algorithm "
405 "needs shell/flexible-constraint support");
407 /* With energy minimization, shells and flexible constraints are
408 * automatically minimized when treated like normal DOFS.
410 if (shellfc != nullptr)
416 if (DOMAINDECOMP(cr))
418 dd_init_local_state(cr->dd, state_global, &ems->s);
420 /* Distribute the charge groups over the nodes from the master node */
421 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
422 imdSession, pull_work, &ems->s, &ems->f, mdAtoms, top, fr, vsite,
423 constr, nrnb, nullptr, FALSE);
424 dd_store_state(cr->dd, &ems->s);
428 state_change_natoms(state_global, state_global->natoms);
429 /* Just copy the state */
430 ems->s = *state_global;
431 state_change_natoms(&ems->s, ems->s.natoms);
433 mdAlgorithmsSetupAtomData(cr, ir, *top_global, top, fr, &ems->f, mdAtoms, constr, vsite,
434 shellfc ? *shellfc : nullptr);
437 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
441 // TODO how should this cross-module support dependency be managed?
442 if (ir->eConstrAlg == econtSHAKE && gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
444 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
445 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
448 if (!ir->bContinuation)
450 /* Constrain the starting coordinates */
451 bool needsLogging = true;
452 bool computeEnergy = true;
453 bool computeVirial = false;
455 constr->apply(needsLogging, computeEnergy, -1, 0, 1.0, ems->s.x.arrayRefWithPadding(),
456 ems->s.x.arrayRefWithPadding(), ArrayRef<RVec>(), ems->s.box,
457 ems->s.lambda[efptFEP], &dvdl_constr, gmx::ArrayRefWithPadding<RVec>(),
458 computeVirial, nullptr, gmx::ConstraintVariable::Positions);
464 *gstat = global_stat_init(ir);
471 calc_shifts(ems->s.box, fr->shift_vec);
474 //! Finalize the minimization
475 static void finish_em(const t_commrec* cr,
477 gmx_walltime_accounting_t walltime_accounting,
478 gmx_wallcycle_t wcycle)
480 if (!thisRankHasDuty(cr, DUTY_PME))
482 /* Tell the PME only node to finish */
483 gmx_pme_send_finish(cr);
488 em_time_end(walltime_accounting, wcycle);
491 //! Swap two different EM states during minimization
492 static void swap_em_state(em_state_t** ems1, em_state_t** ems2)
501 //! Save the EM trajectory
502 static void write_em_traj(FILE* fplog,
508 const gmx_mtop_t* top_global,
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, top_global->natoms, step,
533 static_cast<double>(step), &state->s, state_global,
534 observablesHistory, state->f);
536 if (confout != nullptr)
538 if (DOMAINDECOMP(cr))
540 /* If bX=true, x was collected to state_global in the call above */
543 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
544 dd_collect_vec(cr->dd, &state->s, state->s.x, globalXRef);
549 /* Copy the local state pointer */
550 state_global = &state->s;
555 if (ir->pbcType != PbcType::No && !ir->bPeriodicMols && DOMAINDECOMP(cr))
557 /* Make molecules whole only for confout writing */
558 do_pbc_mtop(ir->pbcType, state->s.box, top_global, state_global->x.rvec_array());
561 write_sto_conf_mtop(confout, *top_global->name, top_global,
562 state_global->x.rvec_array(), nullptr, ir->pbcType, state->s.box);
567 //! \brief Do one minimization step
569 // \returns true when the step succeeded, false when a constraint error occurred
570 static bool do_em_step(const t_commrec* cr,
575 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 /* OpenMP does not supported unsigned loop variables */
664 #pragma omp for schedule(static) nowait
665 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
667 s2->cg_gl[i] = s1->cg_gl[i];
672 if (DOMAINDECOMP(cr))
674 s2->ddp_count = s1->ddp_count;
675 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
681 validStep = constr->apply(
682 TRUE, TRUE, count, 0, 1.0, s1->x.arrayRefWithPadding(), s2->x.arrayRefWithPadding(),
683 ArrayRef<RVec>(), s2->box, s2->lambda[efptBONDED], &dvdl_constr,
684 gmx::ArrayRefWithPadding<RVec>(), false, 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)
701 "The coordinates could not be constrained. Minimizer '%s' can not handle "
702 "constraint failures, use minimizer '%s' before using '%s'.",
703 EI(ir->eI), EI(eiSteep), EI(ir->eI));
710 //! Prepare EM for using domain decomposition parallellization
711 static void em_dd_partition_system(FILE* fplog,
712 const gmx::MDLogger& mdlog,
715 const gmx_mtop_t* top_global,
717 gmx::ImdSession* imdSession,
721 gmx::MDAtoms* mdAtoms,
723 VirtualSitesHandler* vsite,
724 gmx::Constraints* constr,
726 gmx_wallcycle_t wcycle)
728 /* Repartition the domain decomposition */
729 dd_partition_system(fplog, mdlog, step, cr, FALSE, 1, nullptr, *top_global, ir, imdSession, pull_work,
730 &ems->s, &ems->f, mdAtoms, top, fr, vsite, constr, nrnb, wcycle, FALSE);
731 dd_store_state(cr->dd, &ems->s);
737 /*! \brief Class to handle the work of setting and doing an energy evaluation.
739 * This class is a mere aggregate of parameters to pass to evaluate an
740 * energy, so that future changes to names and types of them consume
741 * less time when refactoring other code.
743 * Aggregate initialization is used, for which the chief risk is that
744 * if a member is added at the end and not all initializer lists are
745 * updated, then the member will be value initialized, which will
746 * typically mean initialization to zero.
748 * Use a braced initializer list to construct one of these. */
749 class EnergyEvaluator
752 /*! \brief Evaluates an energy on the state in \c ems.
754 * \todo In practice, the same objects mu_tot, vir, and pres
755 * are always passed to this function, so we would rather have
756 * them as data members. However, their C-array types are
757 * unsuited for aggregate initialization. When the types
758 * improve, the call signature of this method can be reduced.
760 void run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst);
761 //! Handles logging (deprecated).
764 const gmx::MDLogger& mdlog;
765 //! Handles communication.
767 //! Coordinates multi-simulations.
768 const gmx_multisim_t* ms;
769 //! Holds the simulation topology.
770 const gmx_mtop_t* top_global;
771 //! Holds the domain topology.
773 //! User input options.
774 t_inputrec* inputrec;
775 //! The Interactive Molecular Dynamics session.
776 gmx::ImdSession* imdSession;
777 //! The pull work object.
779 //! Manages flop accounting.
781 //! Manages wall cycle accounting.
782 gmx_wallcycle_t wcycle;
783 //! Coordinates global reduction.
784 gmx_global_stat_t gstat;
785 //! Handles virtual sites.
786 VirtualSitesHandler* vsite;
787 //! Handles constraints.
788 gmx::Constraints* constr;
789 //! Per-atom data for this domain.
790 gmx::MDAtoms* mdAtoms;
791 //! Handles how to calculate the forces.
793 //! Schedule of force-calculation work each step for this task.
794 MdrunScheduleWorkload* runScheduleWork;
795 //! Stores the computed energies.
796 gmx_enerdata_t* enerd;
799 void EnergyEvaluator::run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst)
803 tensor force_vir, shake_vir, ekin;
807 /* Set the time to the initial time, the time does not change during EM */
808 t = inputrec->init_t;
810 if (bFirst || (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
812 /* This is the first state or an old state used before the last ns */
818 if (inputrec->nstlist > 0)
826 vsite->construct(ems->s.x, 1, {}, ems->s.box);
829 if (DOMAINDECOMP(cr) && bNS)
831 /* Repartition the domain decomposition */
832 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession, pull_work,
833 ems, top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
836 /* Calc force & energy on new trial position */
837 /* do_force always puts the charge groups in the box and shifts again
838 * We do not unshift, so molecules are always whole in congrad.c
840 do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession, pull_work, count, nrnb, wcycle,
841 top, ems->s.box, ems->s.x.arrayRefWithPadding(), &ems->s.hist,
842 ems->f.arrayRefWithPadding(), force_vir, mdAtoms->mdatoms(), enerd, ems->s.lambda, fr,
843 runScheduleWork, vsite, mu_tot, t, nullptr,
844 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY
845 | (bNS ? GMX_FORCE_NS : 0),
846 DDBalanceRegionHandler(cr));
848 /* Clear the unused shake virial and pressure */
849 clear_mat(shake_vir);
852 /* Communicate stuff when parallel */
853 if (PAR(cr) && inputrec->eI != eiNM)
855 wallcycle_start(wcycle, ewcMoveE);
857 global_stat(gstat, cr, enerd, force_vir, shake_vir, inputrec, nullptr, nullptr, nullptr, 1,
858 &terminate, nullptr, FALSE, CGLO_ENERGY | CGLO_PRESSURE | CGLO_CONSTRAINT);
860 wallcycle_stop(wcycle, ewcMoveE);
863 if (fr->dispersionCorrection)
865 /* Calculate long range corrections to pressure and energy */
866 const DispersionCorrection::Correction correction =
867 fr->dispersionCorrection->calculate(ems->s.box, ems->s.lambda[efptVDW]);
869 enerd->term[F_DISPCORR] = correction.energy;
870 enerd->term[F_EPOT] += correction.energy;
871 enerd->term[F_PRES] += correction.pressure;
872 enerd->term[F_DVDL] += correction.dvdl;
876 enerd->term[F_DISPCORR] = 0;
879 ems->epot = enerd->term[F_EPOT];
883 /* Project out the constraint components of the force */
884 bool needsLogging = false;
885 bool computeEnergy = false;
886 bool computeVirial = true;
888 auto f = ems->f.arrayRefWithPadding();
889 constr->apply(needsLogging, computeEnergy, count, 0, 1.0, ems->s.x.arrayRefWithPadding(), f,
890 f.unpaddedArrayRef(), ems->s.box, ems->s.lambda[efptBONDED], &dvdl_constr,
891 gmx::ArrayRefWithPadding<RVec>(), computeVirial, shake_vir,
892 gmx::ConstraintVariable::ForceDispl);
893 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
894 m_add(force_vir, shake_vir, vir);
898 copy_mat(force_vir, vir);
902 enerd->term[F_PRES] = calc_pres(fr->pbcType, inputrec->nwall, ems->s.box, ekin, vir, pres);
904 sum_dhdl(enerd, ems->s.lambda, *inputrec->fepvals);
906 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
908 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
914 //! Parallel utility summing energies and forces
915 static double reorder_partsum(const t_commrec* cr,
917 const gmx_mtop_t* top_global,
923 fprintf(debug, "Doing reorder_partsum\n");
926 const rvec* fm = s_min->f.rvec_array();
927 const rvec* fb = s_b->f.rvec_array();
929 /* Collect fm in a global vector fmg.
930 * This conflicts with the spirit of domain decomposition,
931 * but to fully optimize this a much more complicated algorithm is required.
933 const int natoms = top_global->natoms;
937 gmx::ArrayRef<const int> indicesMin = s_min->s.cg_gl;
939 for (int a : indicesMin)
941 copy_rvec(fm[i], fmg[a]);
944 gmx_sum(top_global->natoms * 3, fmg[0], cr);
946 /* Now we will determine the part of the sum for the cgs in state s_b */
947 gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
952 gmx::ArrayRef<const unsigned char> grpnrFREEZE =
953 top_global->groups.groupNumbers[SimulationAtomGroupType::Freeze];
954 for (int a : indicesB)
956 if (!grpnrFREEZE.empty())
960 for (int m = 0; m < DIM; m++)
962 if (!opts->nFreeze[gf][m])
964 partsum += (fb[i][m] - fmg[a][m]) * fb[i][m];
975 //! Print some stuff, like beta, whatever that means.
976 static real pr_beta(const t_commrec* cr,
979 const gmx_mtop_t* top_global,
985 /* This is just the classical Polak-Ribiere calculation of beta;
986 * it looks a bit complicated since we take freeze groups into account,
987 * and might have to sum it in parallel runs.
990 if (!DOMAINDECOMP(cr)
991 || (s_min->s.ddp_count == cr->dd->ddp_count && s_b->s.ddp_count == cr->dd->ddp_count))
993 const rvec* fm = s_min->f.rvec_array();
994 const rvec* fb = s_b->f.rvec_array();
997 /* This part of code can be incorrect with DD,
998 * since the atom ordering in s_b and s_min might differ.
1000 for (int i = 0; i < mdatoms->homenr; i++)
1002 if (mdatoms->cFREEZE)
1004 gf = mdatoms->cFREEZE[i];
1006 for (int m = 0; m < DIM; m++)
1008 if (!opts->nFreeze[gf][m])
1010 sum += (fb[i][m] - fm[i][m]) * fb[i][m];
1017 /* We need to reorder cgs while summing */
1018 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1022 gmx_sumd(1, &sum, cr);
1025 return sum / gmx::square(s_min->fnorm);
1031 void LegacySimulator::do_cg()
1033 const char* CG = "Polak-Ribiere Conjugate Gradients";
1035 gmx_localtop_t top(top_global->ffparams);
1036 gmx_global_stat_t gstat;
1037 double tmp, minstep;
1039 real a, b, c, beta = 0.0;
1042 gmx_bool converged, foundlower;
1043 rvec mu_tot = { 0 };
1044 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1046 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1047 int m, step, nminstep;
1048 auto mdatoms = mdAtoms->mdatoms();
1053 "Note that activating conjugate gradient energy minimization via the "
1054 "integrator .mdp option and the command gmx mdrun may "
1055 "be available in a different form in a future version of GROMACS, "
1056 "e.g. gmx minimize and an .mdp option.");
1062 // In CG, the state is extended with a search direction
1063 state_global->flags |= (1 << estCGP);
1065 // Ensure the extra per-atom state array gets allocated
1066 state_change_natoms(state_global, state_global->natoms);
1068 // Initialize the search direction to zero
1069 for (RVec& cg_p : state_global->cg_p)
1075 /* Create 4 states on the stack and extract pointers that we will swap */
1076 em_state_t s0{}, s1{}, s2{}, s3{};
1077 em_state_t* s_min = &s0;
1078 em_state_t* s_a = &s1;
1079 em_state_t* s_b = &s2;
1080 em_state_t* s_c = &s3;
1082 /* Init em and store the local state in s_min */
1083 init_em(fplog, mdlog, CG, cr, inputrec, imdSession, pull_work, state_global, top_global, s_min,
1084 &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, nullptr);
1085 const bool simulationsShareState = false;
1086 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
1087 mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1088 StartingBehavior::NewSimulation, simulationsShareState, ms);
1089 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
1090 false, StartingBehavior::NewSimulation, mdModulesNotifier);
1092 /* Print to log file */
1093 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1095 /* Max number of steps */
1096 number_steps = inputrec->nsteps;
1100 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1104 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1107 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
1108 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1109 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
1110 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1111 /* do_force always puts the charge groups in the box and shifts again
1112 * We do not unshift, so molecules are always whole in congrad.c
1114 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1118 /* Copy stuff to the energy bin for easy printing etc. */
1119 matrix nullBox = {};
1120 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1121 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1122 nullptr, vir, pres, nullptr, mu_tot, constr);
1124 EnergyOutput::printHeader(fplog, step, step);
1125 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step,
1126 step, &fr->listedForces->fcdata(), nullptr);
1129 /* Estimate/guess the initial stepsize */
1130 stepsize = inputrec->em_stepsize / s_min->fnorm;
1134 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1135 fprintf(stderr, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1136 fprintf(stderr, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1137 fprintf(stderr, "\n");
1138 /* and copy to the log file too... */
1139 fprintf(fplog, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1140 fprintf(fplog, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1141 fprintf(fplog, "\n");
1143 /* Start the loop over CG steps.
1144 * Each successful step is counted, and we continue until
1145 * we either converge or reach the max number of steps.
1148 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1151 /* start taking steps in a new direction
1152 * First time we enter the routine, beta=0, and the direction is
1153 * simply the negative gradient.
1156 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1157 rvec* pm = s_min->s.cg_p.rvec_array();
1158 const rvec* sfm = s_min->f.rvec_array();
1161 for (int i = 0; i < mdatoms->homenr; i++)
1163 if (mdatoms->cFREEZE)
1165 gf = mdatoms->cFREEZE[i];
1167 for (m = 0; m < DIM; m++)
1169 if (!inputrec->opts.nFreeze[gf][m])
1171 pm[i][m] = sfm[i][m] + beta * pm[i][m];
1172 gpa -= pm[i][m] * sfm[i][m];
1173 /* f is negative gradient, thus the sign */
1182 /* Sum the gradient along the line across CPUs */
1185 gmx_sumd(1, &gpa, cr);
1188 /* Calculate the norm of the search vector */
1189 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1191 /* Just in case stepsize reaches zero due to numerical precision... */
1194 stepsize = inputrec->em_stepsize / pnorm;
1198 * Double check the value of the derivative in the search direction.
1199 * If it is positive it must be due to the old information in the
1200 * CG formula, so just remove that and start over with beta=0.
1201 * This corresponds to a steepest descent step.
1206 step--; /* Don't count this step since we are restarting */
1207 continue; /* Go back to the beginning of the big for-loop */
1210 /* Calculate minimum allowed stepsize, before the average (norm)
1211 * relative change in coordinate is smaller than precision
1214 auto s_min_x = makeArrayRef(s_min->s.x);
1215 for (int i = 0; i < mdatoms->homenr; i++)
1217 for (m = 0; m < DIM; m++)
1219 tmp = fabs(s_min_x[i][m]);
1224 tmp = pm[i][m] / tmp;
1225 minstep += tmp * tmp;
1228 /* Add up from all CPUs */
1231 gmx_sumd(1, &minstep, cr);
1234 minstep = GMX_REAL_EPS / sqrt(minstep / (3 * top_global->natoms));
1236 if (stepsize < minstep)
1242 /* Write coordinates if necessary */
1243 do_x = do_per_step(step, inputrec->nstxout);
1244 do_f = do_per_step(step, inputrec->nstfout);
1246 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, step, s_min,
1247 state_global, observablesHistory);
1249 /* Take a step downhill.
1250 * In theory, we should minimize the function along this direction.
1251 * That is quite possible, but it turns out to take 5-10 function evaluations
1252 * for each line. However, we dont really need to find the exact minimum -
1253 * it is much better to start a new CG step in a modified direction as soon
1254 * as we are close to it. This will save a lot of energy evaluations.
1256 * In practice, we just try to take a single step.
1257 * If it worked (i.e. lowered the energy), we increase the stepsize but
1258 * the continue straight to the next CG step without trying to find any minimum.
1259 * If it didn't work (higher energy), there must be a minimum somewhere between
1260 * the old position and the new one.
1262 * Due to the finite numerical accuracy, it turns out that it is a good idea
1263 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1264 * This leads to lower final energies in the tests I've done. / Erik
1266 s_a->epot = s_min->epot;
1268 c = a + stepsize; /* reference position along line is zero */
1270 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1272 em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec, imdSession,
1273 pull_work, s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
1276 /* Take a trial step (new coords in s_c) */
1277 do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c, constr, -1);
1280 /* Calculate energy for the trial step */
1281 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1283 /* Calc derivative along line */
1284 const rvec* pc = s_c->s.cg_p.rvec_array();
1285 const rvec* sfc = s_c->f.rvec_array();
1287 for (int i = 0; i < mdatoms->homenr; i++)
1289 for (m = 0; m < DIM; m++)
1291 gpc -= pc[i][m] * sfc[i][m]; /* f is negative gradient, thus the sign */
1294 /* Sum the gradient along the line across CPUs */
1297 gmx_sumd(1, &gpc, cr);
1300 /* This is the max amount of increase in energy we tolerate */
1301 tmp = std::sqrt(GMX_REAL_EPS) * fabs(s_a->epot);
1303 /* Accept the step if the energy is lower, or if it is not significantly higher
1304 * and the line derivative is still negative.
1306 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1309 /* Great, we found a better energy. Increase step for next iteration
1310 * if we are still going down, decrease it otherwise
1314 stepsize *= 1.618034; /* The golden section */
1318 stepsize *= 0.618034; /* 1/golden section */
1323 /* New energy is the same or higher. We will have to do some work
1324 * to find a smaller value in the interval. Take smaller step next time!
1327 stepsize *= 0.618034;
1331 /* OK, if we didn't find a lower value we will have to locate one now - there must
1332 * be one in the interval [a=0,c].
1333 * The same thing is valid here, though: Don't spend dozens of iterations to find
1334 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1335 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1337 * I also have a safeguard for potentially really pathological functions so we never
1338 * take more than 20 steps before we give up ...
1340 * If we already found a lower value we just skip this step and continue to the update.
1349 /* Select a new trial point.
1350 * If the derivatives at points a & c have different sign we interpolate to zero,
1351 * otherwise just do a bisection.
1353 if (gpa < 0 && gpc > 0)
1355 b = a + gpa * (a - c) / (gpc - gpa);
1362 /* safeguard if interpolation close to machine accuracy causes errors:
1363 * never go outside the interval
1365 if (b <= a || b >= c)
1370 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1372 /* Reload the old state */
1373 em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec, imdSession, pull_work,
1374 s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
1377 /* Take a trial step to this new point - new coords in s_b */
1378 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b, constr, -1);
1381 /* Calculate energy for the trial step */
1382 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1384 /* p does not change within a step, but since the domain decomposition
1385 * might change, we have to use cg_p of s_b here.
1387 const rvec* pb = s_b->s.cg_p.rvec_array();
1388 const rvec* sfb = s_b->f.rvec_array();
1390 for (int i = 0; i < mdatoms->homenr; i++)
1392 for (m = 0; m < DIM; m++)
1394 gpb -= pb[i][m] * sfb[i][m]; /* f is negative gradient, thus the sign */
1397 /* Sum the gradient along the line across CPUs */
1400 gmx_sumd(1, &gpb, cr);
1405 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a->epot, s_b->epot,
1409 epot_repl = s_b->epot;
1411 /* Keep one of the intervals based on the value of the derivative at the new point */
1414 /* Replace c endpoint with b */
1415 swap_em_state(&s_b, &s_c);
1421 /* Replace a endpoint with b */
1422 swap_em_state(&s_b, &s_a);
1428 * Stop search as soon as we find a value smaller than the endpoints.
1429 * Never run more than 20 steps, no matter what.
1432 } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) && (nminstep < 20));
1434 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot) * GMX_REAL_EPS || nminstep >= 20)
1436 /* OK. We couldn't find a significantly lower energy.
1437 * If beta==0 this was steepest descent, and then we give up.
1438 * If not, set beta=0 and restart with steepest descent before quitting.
1448 /* Reset memory before giving up */
1454 /* Select min energy state of A & C, put the best in B.
1456 if (s_c->epot < s_a->epot)
1460 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c->epot,
1463 swap_em_state(&s_b, &s_c);
1470 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a->epot,
1473 swap_em_state(&s_b, &s_a);
1481 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n", s_c->epot);
1483 swap_em_state(&s_b, &s_c);
1487 /* new search direction */
1488 /* beta = 0 means forget all memory and restart with steepest descents. */
1489 if (nstcg && ((step % nstcg) == 0))
1495 /* s_min->fnorm cannot be zero, because then we would have converged
1499 /* Polak-Ribiere update.
1500 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1502 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1504 /* Limit beta to prevent oscillations */
1505 if (fabs(beta) > 5.0)
1511 /* update positions */
1512 swap_em_state(&s_min, &s_b);
1515 /* Print it if necessary */
1518 if (mdrunOptions.verbose)
1520 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1521 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step,
1522 s_min->epot, s_min->fnorm / sqrtNumAtoms, s_min->fmax, s_min->a_fmax + 1);
1525 /* Store the new (lower) energies */
1526 matrix nullBox = {};
1527 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1528 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1529 nullptr, vir, pres, nullptr, mu_tot, constr);
1531 do_log = do_per_step(step, inputrec->nstlog);
1532 do_ene = do_per_step(step, inputrec->nstenergy);
1534 imdSession->fillEnergyRecord(step, TRUE);
1538 EnergyOutput::printHeader(fplog, step, step);
1540 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1541 do_log ? fplog : nullptr, step, step,
1542 &fr->listedForces->fcdata(), nullptr);
1545 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1546 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0))
1548 imdSession->sendPositionsAndEnergies();
1551 /* Stop when the maximum force lies below tolerance.
1552 * If we have reached machine precision, converged is already set to true.
1554 converged = converged || (s_min->fmax < inputrec->em_tol);
1556 } /* End of the loop */
1560 step--; /* we never took that last step in this case */
1562 if (s_min->fmax > inputrec->em_tol)
1566 warn_step(fplog, inputrec->em_tol, s_min->fmax, step - 1 == number_steps, FALSE);
1573 /* If we printed energy and/or logfile last step (which was the last step)
1574 * we don't have to do it again, but otherwise print the final values.
1578 /* Write final value to log since we didn't do anything the last step */
1579 EnergyOutput::printHeader(fplog, step, step);
1581 if (!do_ene || !do_log)
1583 /* Write final energy file entries */
1584 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1585 !do_log ? fplog : nullptr, step, step,
1586 &fr->listedForces->fcdata(), nullptr);
1590 /* Print some stuff... */
1593 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1597 * For accurate normal mode calculation it is imperative that we
1598 * store the last conformation into the full precision binary trajectory.
1600 * However, we should only do it if we did NOT already write this step
1601 * above (which we did if do_x or do_f was true).
1603 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1604 * in the trajectory with the same step number.
1606 do_x = !do_per_step(step, inputrec->nstxout);
1607 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1609 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec,
1610 step, s_min, state_global, observablesHistory);
1615 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1616 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1617 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1619 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1622 finish_em(cr, outf, walltime_accounting, wcycle);
1624 /* To print the actual number of steps we needed somewhere */
1625 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1629 void LegacySimulator::do_lbfgs()
1631 static const char* LBFGS = "Low-Memory BFGS Minimizer";
1633 gmx_localtop_t top(top_global->ffparams);
1634 gmx_global_stat_t gstat;
1635 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1636 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1637 real * rho, *alpha, *p, *s, **dx, **dg;
1638 real a, b, c, maxdelta, delta;
1640 real dgdx, dgdg, sq, yr, beta;
1642 rvec mu_tot = { 0 };
1643 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1645 int start, end, number_steps;
1646 int i, k, m, n, gf, step;
1648 auto mdatoms = mdAtoms->mdatoms();
1653 "Note that activating L-BFGS energy minimization via the "
1654 "integrator .mdp option and the command gmx mdrun may "
1655 "be available in a different form in a future version of GROMACS, "
1656 "e.g. gmx minimize and an .mdp option.");
1660 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1663 if (nullptr != constr)
1667 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1668 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1671 n = 3 * state_global->natoms;
1672 nmaxcorr = inputrec->nbfgscorr;
1677 snew(rho, nmaxcorr);
1678 snew(alpha, nmaxcorr);
1681 for (i = 0; i < nmaxcorr; i++)
1687 for (i = 0; i < nmaxcorr; i++)
1696 init_em(fplog, mdlog, LBFGS, cr, inputrec, imdSession, pull_work, state_global, top_global,
1697 &ems, &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, nullptr);
1698 const bool simulationsShareState = false;
1699 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
1700 mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1701 StartingBehavior::NewSimulation, simulationsShareState, ms);
1702 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
1703 false, StartingBehavior::NewSimulation, mdModulesNotifier);
1706 end = mdatoms->homenr;
1708 /* We need 4 working states */
1709 em_state_t s0{}, s1{}, s2{}, s3{};
1710 em_state_t* sa = &s0;
1711 em_state_t* sb = &s1;
1712 em_state_t* sc = &s2;
1713 em_state_t* last = &s3;
1714 /* Initialize by copying the state from ems (we could skip x and f here) */
1719 /* Print to log file */
1720 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1722 do_log = do_ene = do_x = do_f = TRUE;
1724 /* Max number of steps */
1725 number_steps = inputrec->nsteps;
1727 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1729 for (i = start; i < end; i++)
1731 if (mdatoms->cFREEZE)
1733 gf = mdatoms->cFREEZE[i];
1735 for (m = 0; m < DIM; m++)
1737 frozen[3 * i + m] = (inputrec->opts.nFreeze[gf][m] != 0);
1742 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1746 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1751 vsite->construct(state_global->x, 1, {}, state_global->box);
1754 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1755 /* do_force always puts the charge groups in the box and shifts again
1756 * We do not unshift, so molecules are always whole
1759 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
1760 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1761 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
1762 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1766 /* Copy stuff to the energy bin for easy printing etc. */
1767 matrix nullBox = {};
1768 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1769 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1770 nullptr, vir, pres, nullptr, mu_tot, constr);
1772 EnergyOutput::printHeader(fplog, step, step);
1773 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step,
1774 step, &fr->listedForces->fcdata(), nullptr);
1777 /* Set the initial step.
1778 * since it will be multiplied by the non-normalized search direction
1779 * vector (force vector the first time), we scale it by the
1780 * norm of the force.
1785 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1786 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1787 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1788 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
1789 fprintf(stderr, "\n");
1790 /* and copy to the log file too... */
1791 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1792 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1793 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
1794 fprintf(fplog, "\n");
1797 // Point is an index to the memory of search directions, where 0 is the first one.
1800 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1801 real* fInit = static_cast<real*>(ems.f.rvec_array()[0]);
1802 for (i = 0; i < n; i++)
1806 dx[point][i] = fInit[i]; /* Initial search direction */
1814 // Stepsize will be modified during the search, and actually it is not critical
1815 // (the main efficiency in the algorithm comes from changing directions), but
1816 // we still need an initial value, so estimate it as the inverse of the norm
1817 // so we take small steps where the potential fluctuates a lot.
1818 stepsize = 1.0 / ems.fnorm;
1820 /* Start the loop over BFGS steps.
1821 * Each successful step is counted, and we continue until
1822 * we either converge or reach the max number of steps.
1827 /* Set the gradient from the force */
1829 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1832 /* Write coordinates if necessary */
1833 do_x = do_per_step(step, inputrec->nstxout);
1834 do_f = do_per_step(step, inputrec->nstfout);
1839 mdof_flags |= MDOF_X;
1844 mdof_flags |= MDOF_F;
1849 mdof_flags |= MDOF_IMD;
1852 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global->natoms, step,
1853 static_cast<real>(step), &ems.s, state_global,
1854 observablesHistory, ems.f);
1856 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1858 /* make s a pointer to current search direction - point=0 first time we get here */
1861 real* xx = static_cast<real*>(ems.s.x.rvec_array()[0]);
1862 real* ff = static_cast<real*>(ems.f.rvec_array()[0]);
1864 // calculate line gradient in position A
1865 for (gpa = 0, i = 0; i < n; i++)
1867 gpa -= s[i] * ff[i];
1870 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1871 * relative change in coordinate is smaller than precision
1873 for (minstep = 0, i = 0; i < n; i++)
1881 minstep += tmp * tmp;
1883 minstep = GMX_REAL_EPS / sqrt(minstep / n);
1885 if (stepsize < minstep)
1891 // Before taking any steps along the line, store the old position
1893 real* lastx = static_cast<real*>(last->s.x.data()[0]);
1894 real* lastf = static_cast<real*>(last->f.data()[0]);
1899 /* Take a step downhill.
1900 * In theory, we should find the actual minimum of the function in this
1901 * direction, somewhere along the line.
1902 * That is quite possible, but it turns out to take 5-10 function evaluations
1903 * for each line. However, we dont really need to find the exact minimum -
1904 * it is much better to start a new BFGS step in a modified direction as soon
1905 * as we are close to it. This will save a lot of energy evaluations.
1907 * In practice, we just try to take a single step.
1908 * If it worked (i.e. lowered the energy), we increase the stepsize but
1909 * continue straight to the next BFGS step without trying to find any minimum,
1910 * i.e. we change the search direction too. If the line was smooth, it is
1911 * likely we are in a smooth region, and then it makes sense to take longer
1912 * steps in the modified search direction too.
1914 * If it didn't work (higher energy), there must be a minimum somewhere between
1915 * the old position and the new one. Then we need to start by finding a lower
1916 * value before we change search direction. Since the energy was apparently
1917 * quite rough, we need to decrease the step size.
1919 * Due to the finite numerical accuracy, it turns out that it is a good idea
1920 * to accept a SMALL increase in energy, if the derivative is still downhill.
1921 * This leads to lower final energies in the tests I've done. / Erik
1924 // State "A" is the first position along the line.
1925 // reference position along line is initially zero
1928 // Check stepsize first. We do not allow displacements
1929 // larger than emstep.
1933 // Pick a new position C by adding stepsize to A.
1936 // Calculate what the largest change in any individual coordinate
1937 // would be (translation along line * gradient along line)
1939 for (i = 0; i < n; i++)
1942 if (delta > maxdelta)
1947 // If any displacement is larger than the stepsize limit, reduce the step
1948 if (maxdelta > inputrec->em_stepsize)
1952 } while (maxdelta > inputrec->em_stepsize);
1954 // Take a trial step and move the coordinate array xc[] to position C
1955 real* xc = static_cast<real*>(sc->s.x.rvec_array()[0]);
1956 for (i = 0; i < n; i++)
1958 xc[i] = lastx[i] + c * s[i];
1962 // Calculate energy for the trial step in position C
1963 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
1965 // Calc line gradient in position C
1966 real* fc = static_cast<real*>(sc->f.rvec_array()[0]);
1967 for (gpc = 0, i = 0; i < n; i++)
1969 gpc -= s[i] * fc[i]; /* f is negative gradient, thus the sign */
1971 /* Sum the gradient along the line across CPUs */
1974 gmx_sumd(1, &gpc, cr);
1977 // This is the max amount of increase in energy we tolerate.
1978 // By allowing VERY small changes (close to numerical precision) we
1979 // frequently find even better (lower) final energies.
1980 tmp = std::sqrt(GMX_REAL_EPS) * fabs(sa->epot);
1982 // Accept the step if the energy is lower in the new position C (compared to A),
1983 // or if it is not significantly higher and the line derivative is still negative.
1984 foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
1985 // If true, great, we found a better energy. We no longer try to alter the
1986 // stepsize, but simply accept this new better position. The we select a new
1987 // search direction instead, which will be much more efficient than continuing
1988 // to take smaller steps along a line. Set fnorm based on the new C position,
1989 // which will be used to update the stepsize to 1/fnorm further down.
1991 // If false, the energy is NOT lower in point C, i.e. it will be the same
1992 // or higher than in point A. In this case it is pointless to move to point C,
1993 // so we will have to do more iterations along the same line to find a smaller
1994 // value in the interval [A=0.0,C].
1995 // Here, A is still 0.0, but that will change when we do a search in the interval
1996 // [0.0,C] below. That search we will do by interpolation or bisection rather
1997 // than with the stepsize, so no need to modify it. For the next search direction
1998 // it will be reset to 1/fnorm anyway.
2002 // OK, if we didn't find a lower value we will have to locate one now - there must
2003 // be one in the interval [a,c].
2004 // The same thing is valid here, though: Don't spend dozens of iterations to find
2005 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2006 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2007 // I also have a safeguard for potentially really pathological functions so we never
2008 // take more than 20 steps before we give up.
2009 // If we already found a lower value we just skip this step and continue to the update.
2014 // Select a new trial point B in the interval [A,C].
2015 // If the derivatives at points a & c have different sign we interpolate to zero,
2016 // otherwise just do a bisection since there might be multiple minima/maxima
2017 // inside the interval.
2018 if (gpa < 0 && gpc > 0)
2020 b = a + gpa * (a - c) / (gpc - gpa);
2027 /* safeguard if interpolation close to machine accuracy causes errors:
2028 * never go outside the interval
2030 if (b <= a || b >= c)
2035 // Take a trial step to point B
2036 real* xb = static_cast<real*>(sb->s.x.rvec_array()[0]);
2037 for (i = 0; i < n; i++)
2039 xb[i] = lastx[i] + b * s[i];
2043 // Calculate energy for the trial step in point B
2044 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2047 // Calculate gradient in point B
2048 real* fb = static_cast<real*>(sb->f.rvec_array()[0]);
2049 for (gpb = 0, i = 0; i < n; i++)
2051 gpb -= s[i] * fb[i]; /* f is negative gradient, thus the sign */
2053 /* Sum the gradient along the line across CPUs */
2056 gmx_sumd(1, &gpb, cr);
2059 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2060 // at the new point B, and rename the endpoints of this new interval A and C.
2063 /* Replace c endpoint with b */
2065 /* copy state b to c */
2070 /* Replace a endpoint with b */
2072 /* copy state b to a */
2077 * Stop search as soon as we find a value smaller than the endpoints,
2078 * or if the tolerance is below machine precision.
2079 * Never run more than 20 steps, no matter what.
2082 } while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2084 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2086 /* OK. We couldn't find a significantly lower energy.
2087 * If ncorr==0 this was steepest descent, and then we give up.
2088 * If not, reset memory to restart as steepest descent before quitting.
2100 /* Search in gradient direction */
2101 for (i = 0; i < n; i++)
2103 dx[point][i] = ff[i];
2105 /* Reset stepsize */
2106 stepsize = 1.0 / fnorm;
2111 /* Select min energy state of A & C, put the best in xx/ff/Epot
2113 if (sc->epot < sa->epot)
2134 /* Update the memory information, and calculate a new
2135 * approximation of the inverse hessian
2138 /* Have new data in Epot, xx, ff */
2139 if (ncorr < nmaxcorr)
2144 for (i = 0; i < n; i++)
2146 dg[point][i] = lastf[i] - ff[i];
2147 dx[point][i] *= step_taken;
2152 for (i = 0; i < n; i++)
2154 dgdg += dg[point][i] * dg[point][i];
2155 dgdx += dg[point][i] * dx[point][i];
2160 rho[point] = 1.0 / dgdx;
2163 if (point >= nmaxcorr)
2169 for (i = 0; i < n; i++)
2176 /* Recursive update. First go back over the memory points */
2177 for (k = 0; k < ncorr; k++)
2186 for (i = 0; i < n; i++)
2188 sq += dx[cp][i] * p[i];
2191 alpha[cp] = rho[cp] * sq;
2193 for (i = 0; i < n; i++)
2195 p[i] -= alpha[cp] * dg[cp][i];
2199 for (i = 0; i < n; i++)
2204 /* And then go forward again */
2205 for (k = 0; k < ncorr; k++)
2208 for (i = 0; i < n; i++)
2210 yr += p[i] * dg[cp][i];
2213 beta = rho[cp] * yr;
2214 beta = alpha[cp] - beta;
2216 for (i = 0; i < n; i++)
2218 p[i] += beta * dx[cp][i];
2228 for (i = 0; i < n; i++)
2232 dx[point][i] = p[i];
2240 /* Print it if necessary */
2243 if (mdrunOptions.verbose)
2245 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2246 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step,
2247 ems.epot, ems.fnorm / sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2250 /* Store the new (lower) energies */
2251 matrix nullBox = {};
2252 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
2253 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
2254 nullptr, vir, pres, nullptr, mu_tot, constr);
2256 do_log = do_per_step(step, inputrec->nstlog);
2257 do_ene = do_per_step(step, inputrec->nstenergy);
2259 imdSession->fillEnergyRecord(step, TRUE);
2263 EnergyOutput::printHeader(fplog, step, step);
2265 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2266 do_log ? fplog : nullptr, step, step,
2267 &fr->listedForces->fcdata(), nullptr);
2270 /* Send x and E to IMD client, if bIMD is TRUE. */
2271 if (imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0) && MASTER(cr))
2273 imdSession->sendPositionsAndEnergies();
2276 // Reset stepsize in we are doing more iterations
2279 /* Stop when the maximum force lies below tolerance.
2280 * If we have reached machine precision, converged is already set to true.
2282 converged = converged || (ems.fmax < inputrec->em_tol);
2284 } /* End of the loop */
2288 step--; /* we never took that last step in this case */
2290 if (ems.fmax > inputrec->em_tol)
2294 warn_step(fplog, inputrec->em_tol, ems.fmax, step - 1 == number_steps, FALSE);
2299 /* If we printed energy and/or logfile last step (which was the last step)
2300 * we don't have to do it again, but otherwise print the final values.
2302 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2304 EnergyOutput::printHeader(fplog, step, step);
2306 if (!do_ene || !do_log) /* Write final energy file entries */
2308 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2309 !do_log ? fplog : nullptr, step, step,
2310 &fr->listedForces->fcdata(), nullptr);
2313 /* Print some stuff... */
2316 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2320 * For accurate normal mode calculation it is imperative that we
2321 * store the last conformation into the full precision binary trajectory.
2323 * However, we should only do it if we did NOT already write this step
2324 * above (which we did if do_x or do_f was true).
2326 do_x = !do_per_step(step, inputrec->nstxout);
2327 do_f = !do_per_step(step, inputrec->nstfout);
2328 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec,
2329 step, &ems, state_global, observablesHistory);
2333 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2334 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2335 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2337 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2340 finish_em(cr, outf, walltime_accounting, wcycle);
2342 /* To print the actual number of steps we needed somewhere */
2343 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2346 void LegacySimulator::do_steep()
2348 const char* SD = "Steepest Descents";
2349 gmx_localtop_t top(top_global->ffparams);
2350 gmx_global_stat_t gstat;
2353 gmx_bool bDone, bAbort, do_x, do_f;
2355 rvec mu_tot = { 0 };
2358 int steps_accepted = 0;
2359 auto mdatoms = mdAtoms->mdatoms();
2364 "Note that activating steepest-descent energy minimization via the "
2365 "integrator .mdp option and the command gmx mdrun may "
2366 "be available in a different form in a future version of GROMACS, "
2367 "e.g. gmx minimize and an .mdp option.");
2369 /* Create 2 states on the stack and extract pointers that we will swap */
2370 em_state_t s0{}, s1{};
2371 em_state_t* s_min = &s0;
2372 em_state_t* s_try = &s1;
2374 /* Init em and store the local state in s_try */
2375 init_em(fplog, mdlog, SD, cr, inputrec, imdSession, pull_work, state_global, top_global, s_try,
2376 &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, nullptr);
2377 const bool simulationsShareState = false;
2378 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
2379 mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2380 StartingBehavior::NewSimulation, simulationsShareState, ms);
2381 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
2382 false, StartingBehavior::NewSimulation, mdModulesNotifier);
2384 /* Print to log file */
2385 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2387 /* Set variables for stepsize (in nm). This is the largest
2388 * step that we are going to make in any direction.
2390 ustep = inputrec->em_stepsize;
2393 /* Max number of steps */
2394 nsteps = inputrec->nsteps;
2398 /* Print to the screen */
2399 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2403 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2405 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2406 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2407 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2409 /**** HERE STARTS THE LOOP ****
2410 * count is the counter for the number of steps
2411 * bDone will be TRUE when the minimization has converged
2412 * bAbort will be TRUE when nsteps steps have been performed or when
2413 * the stepsize becomes smaller than is reasonable for machine precision
2418 while (!bDone && !bAbort)
2420 bAbort = (nsteps >= 0) && (count == nsteps);
2422 /* set new coordinates, except for first step */
2423 bool validStep = true;
2426 validStep = do_em_step(cr, inputrec, mdatoms, s_min, stepsize, &s_min->f, s_try, constr, count);
2431 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2435 // Signal constraint error during stepping with energy=inf
2436 s_try->epot = std::numeric_limits<real>::infinity();
2441 EnergyOutput::printHeader(fplog, count, count);
2446 s_min->epot = s_try->epot;
2449 /* Print it if necessary */
2452 if (mdrunOptions.verbose)
2454 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2455 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax + 1,
2456 ((count == 0) || (s_try->epot < s_min->epot)) ? '\n' : '\r');
2460 if ((count == 0) || (s_try->epot < s_min->epot))
2462 /* Store the new (lower) energies */
2463 matrix nullBox = {};
2464 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(count), mdatoms->tmass,
2465 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
2466 nullptr, vir, pres, nullptr, mu_tot, constr);
2468 imdSession->fillEnergyRecord(count, TRUE);
2470 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2471 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2472 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, do_dr, do_or, fplog,
2473 count, count, &fr->listedForces->fcdata(), nullptr);
2478 /* Now if the new energy is smaller than the previous...
2479 * or if this is the first step!
2480 * or if we did random steps!
2483 if ((count == 0) || (s_try->epot < s_min->epot))
2487 /* Test whether the convergence criterion is met... */
2488 bDone = (s_try->fmax < inputrec->em_tol);
2490 /* Copy the arrays for force, positions and energy */
2491 /* The 'Min' array always holds the coords and forces of the minimal
2493 swap_em_state(&s_min, &s_try);
2499 /* Write to trn, if necessary */
2500 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2501 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2502 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, count, s_min,
2503 state_global, observablesHistory);
2507 /* If energy is not smaller make the step smaller... */
2510 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2512 /* Reload the old state */
2513 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession,
2514 pull_work, s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
2518 // If the force is very small after finishing minimization,
2519 // we risk dividing by zero when calculating the step size.
2520 // So we check first if the minimization has stopped before
2521 // trying to obtain a new step size.
2524 /* Determine new step */
2525 stepsize = ustep / s_min->fmax;
2528 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2530 if (count == nsteps || ustep < 1e-12)
2532 if (count == nsteps || ustep < 1e-6)
2537 warn_step(fplog, inputrec->em_tol, s_min->fmax, count == nsteps, constr != nullptr);
2542 /* Send IMD energies and positions, if bIMD is TRUE. */
2543 if (imdSession->run(count, TRUE, state_global->box,
2544 MASTER(cr) ? state_global->x.rvec_array() : nullptr, 0)
2547 imdSession->sendPositionsAndEnergies();
2551 } /* End of the loop */
2553 /* Print some data... */
2556 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2558 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2559 top_global, inputrec, count, s_min, state_global, observablesHistory);
2563 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2565 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2566 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2569 finish_em(cr, outf, walltime_accounting, wcycle);
2571 /* To print the actual number of steps we needed somewhere */
2572 inputrec->nsteps = count;
2574 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2577 void LegacySimulator::do_nm()
2579 const char* NM = "Normal Mode Analysis";
2581 gmx_localtop_t top(top_global->ffparams);
2582 gmx_global_stat_t gstat;
2584 rvec mu_tot = { 0 };
2586 gmx_bool bSparse; /* use sparse matrix storage format */
2588 gmx_sparsematrix_t* sparse_matrix = nullptr;
2589 real* full_matrix = nullptr;
2591 /* added with respect to mdrun */
2593 real der_range = 10.0 * std::sqrt(GMX_REAL_EPS);
2595 bool bIsMaster = MASTER(cr);
2596 auto mdatoms = mdAtoms->mdatoms();
2601 "Note that activating normal-mode analysis via the integrator "
2602 ".mdp option and the command gmx mdrun may "
2603 "be available in a different form in a future version of GROMACS, "
2604 "e.g. gmx normal-modes.");
2606 if (constr != nullptr)
2610 "Constraints present with Normal Mode Analysis, this combination is not supported");
2613 gmx_shellfc_t* shellfc;
2615 em_state_t state_work{};
2617 /* Init em and store the local state in state_minimum */
2618 init_em(fplog, mdlog, NM, cr, inputrec, imdSession, pull_work, state_global, top_global,
2619 &state_work, &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, &shellfc);
2620 const bool simulationsShareState = false;
2621 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
2622 mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2623 StartingBehavior::NewSimulation, simulationsShareState, ms);
2625 std::vector<int> atom_index = get_atom_index(top_global);
2626 std::vector<gmx::RVec> fneg(atom_index.size(), { 0, 0, 0 });
2627 snew(dfdx, atom_index.size());
2633 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2634 " which MIGHT not be accurate enough for normal mode analysis.\n"
2635 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2636 " are fairly modest even if you recompile in double precision.\n\n");
2640 /* Check if we can/should use sparse storage format.
2642 * Sparse format is only useful when the Hessian itself is sparse, which it
2643 * will be when we use a cutoff.
2644 * For small systems (n<1000) it is easier to always use full matrix format, though.
2646 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2648 GMX_LOG(mdlog.warning)
2649 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2652 else if (atom_index.size() < 1000)
2654 GMX_LOG(mdlog.warning)
2655 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2661 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2665 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2666 sz = DIM * atom_index.size();
2668 fprintf(stderr, "Allocating Hessian memory...\n\n");
2672 sparse_matrix = gmx_sparsematrix_init(sz);
2673 sparse_matrix->compressed_symmetric = TRUE;
2677 snew(full_matrix, sz * sz);
2680 /* Write start time and temperature */
2681 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2683 /* fudge nr of steps to nr of atoms */
2684 inputrec->nsteps = atom_index.size() * 2;
2688 fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2689 *(top_global->name), inputrec->nsteps);
2692 nnodes = cr->nnodes;
2694 /* Make evaluate_energy do a single node force calculation */
2696 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2697 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2698 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2699 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2700 cr->nnodes = nnodes;
2702 /* if forces are not small, warn user */
2703 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2705 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2706 if (state_work.fmax > 1.0e-3)
2708 GMX_LOG(mdlog.warning)
2710 "The force is probably not small enough to "
2711 "ensure that you are at a minimum.\n"
2712 "Be aware that negative eigenvalues may occur\n"
2713 "when the resulting matrix is diagonalized.");
2716 /***********************************************************
2718 * Loop over all pairs in matrix
2720 * do_force called twice. Once with positive and
2721 * once with negative displacement
2723 ************************************************************/
2725 /* Steps are divided one by one over the nodes */
2727 auto state_work_x = makeArrayRef(state_work.s.x);
2728 auto state_work_f = makeArrayRef(state_work.f);
2729 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
2731 size_t atom = atom_index[aid];
2732 for (size_t d = 0; d < DIM; d++)
2735 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2738 x_min = state_work_x[atom][d];
2740 for (unsigned int dx = 0; (dx < 2); dx++)
2744 state_work_x[atom][d] = x_min - der_range;
2748 state_work_x[atom][d] = x_min + der_range;
2751 /* Make evaluate_energy do a single node force calculation */
2755 /* Now is the time to relax the shells */
2756 relax_shell_flexcon(
2757 fplog, cr, ms, mdrunOptions.verbose, nullptr, step, inputrec, imdSession,
2758 pull_work, bNS, force_flags, &top, constr, enerd, state_work.s.natoms,
2759 state_work.s.x.arrayRefWithPadding(), state_work.s.v.arrayRefWithPadding(),
2760 state_work.s.box, state_work.s.lambda, &state_work.s.hist,
2761 state_work.f.arrayRefWithPadding(), vir, mdatoms, nrnb, wcycle, shellfc,
2762 fr, runScheduleWork, t, mu_tot, vsite, DDBalanceRegionHandler(nullptr));
2768 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid * 2 + dx, FALSE);
2771 cr->nnodes = nnodes;
2775 std::copy(state_work_f.begin(), state_work_f.begin() + atom_index.size(),
2780 /* x is restored to original */
2781 state_work_x[atom][d] = x_min;
2783 for (size_t j = 0; j < atom_index.size(); j++)
2785 for (size_t k = 0; (k < DIM); k++)
2787 dfdx[j][k] = -(state_work_f[atom_index[j]][k] - fneg[j][k]) / (2 * der_range);
2794 # define mpi_type GMX_MPI_REAL
2795 MPI_Send(dfdx[0], atom_index.size() * DIM, mpi_type, MASTER(cr), cr->nodeid,
2796 cr->mpi_comm_mygroup);
2801 for (index node = 0; (node < nnodes && aid + node < ssize(atom_index)); node++)
2807 MPI_Recv(dfdx[0], atom_index.size() * DIM, mpi_type, node, node,
2808 cr->mpi_comm_mygroup, &stat);
2813 row = (aid + node) * DIM + d;
2815 for (size_t j = 0; j < atom_index.size(); j++)
2817 for (size_t k = 0; k < DIM; k++)
2823 if (col >= row && dfdx[j][k] != 0.0)
2825 gmx_sparsematrix_increment_value(sparse_matrix, row, col, dfdx[j][k]);
2830 full_matrix[row * sz + col] = dfdx[j][k];
2837 if (mdrunOptions.verbose && fplog)
2842 /* write progress */
2843 if (bIsMaster && mdrunOptions.verbose)
2845 fprintf(stderr, "\rFinished step %d out of %td",
2846 std::min<int>(atom + nnodes, atom_index.size()), ssize(atom_index));
2853 fprintf(stderr, "\n\nWriting Hessian...\n");
2854 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2857 finish_em(cr, outf, walltime_accounting, wcycle);
2859 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size() * 2);