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,2021, 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/checkpointdata.h"
94 #include "gromacs/mdtypes/commrec.h"
95 #include "gromacs/mdtypes/forcebuffers.h"
96 #include "gromacs/mdtypes/forcerec.h"
97 #include "gromacs/mdtypes/inputrec.h"
98 #include "gromacs/mdtypes/interaction_const.h"
99 #include "gromacs/mdtypes/md_enums.h"
100 #include "gromacs/mdtypes/mdatom.h"
101 #include "gromacs/mdtypes/mdrunoptions.h"
102 #include "gromacs/mdtypes/state.h"
103 #include "gromacs/pbcutil/pbc.h"
104 #include "gromacs/timing/wallcycle.h"
105 #include "gromacs/timing/walltime_accounting.h"
106 #include "gromacs/topology/mtop_util.h"
107 #include "gromacs/topology/topology.h"
108 #include "gromacs/utility/cstringutil.h"
109 #include "gromacs/utility/exceptions.h"
110 #include "gromacs/utility/fatalerror.h"
111 #include "gromacs/utility/logger.h"
112 #include "gromacs/utility/smalloc.h"
114 #include "legacysimulator.h"
118 using gmx::MdrunScheduleWorkload;
120 using gmx::VirtualSitesHandler;
122 //! Utility structure for manipulating states during EM
123 typedef struct em_state
125 //! Copy of the global state
131 //! Norm of the force
139 //! Print the EM starting conditions
140 static void print_em_start(FILE* fplog,
142 gmx_walltime_accounting_t walltime_accounting,
143 gmx_wallcycle* wcycle,
146 walltime_accounting_start_time(walltime_accounting);
147 wallcycle_start(wcycle, WallCycleCounter::Run);
148 print_start(fplog, cr, walltime_accounting, name);
151 //! Stop counting time for EM
152 static void em_time_end(gmx_walltime_accounting_t walltime_accounting, gmx_wallcycle* wcycle)
154 wallcycle_stop(wcycle, WallCycleCounter::Run);
156 walltime_accounting_end_time(walltime_accounting);
159 //! Printing a log file and console header
160 static void sp_header(FILE* out, const char* minimizer, real ftol, int nsteps)
163 fprintf(out, "%s:\n", minimizer);
164 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
165 fprintf(out, " Number of steps = %12d\n", nsteps);
168 //! Print warning message
169 static void warn_step(FILE* fp, real ftol, real fmax, gmx_bool bLastStep, gmx_bool bConstrain)
171 constexpr bool realIsDouble = GMX_DOUBLE;
174 if (!std::isfinite(fmax))
177 "\nEnergy minimization has stopped because the force "
178 "on at least one atom is not finite. This usually means "
179 "atoms are overlapping. Modify the input coordinates to "
180 "remove atom overlap or use soft-core potentials with "
181 "the free energy code to avoid infinite forces.\n%s",
182 !realIsDouble ? "You could also be lucky that switching to double precision "
183 "is sufficient to obtain finite forces.\n"
189 "\nEnergy minimization reached the maximum number "
190 "of steps before the forces reached the requested "
191 "precision Fmax < %g.\n",
197 "\nEnergy minimization has stopped, but the forces have "
198 "not converged to the requested precision Fmax < %g (which "
199 "may not be possible for your system). It stopped "
200 "because the algorithm tried to make a new step whose size "
201 "was too small, or there was no change in the energy since "
202 "last step. Either way, we regard the minimization as "
203 "converged to within the available machine precision, "
204 "given your starting configuration and EM parameters.\n%s%s",
206 !realIsDouble ? "\nDouble precision normally gives you higher accuracy, but "
207 "this is often not needed for preparing to run molecular "
210 bConstrain ? "You might need to increase your constraint accuracy, or turn\n"
211 "off constraints altogether (set constraints = none in mdp file)\n"
215 fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
216 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
219 //! Print message about convergence of the EM
220 static void print_converged(FILE* fp,
226 const em_state_t* ems,
229 char buf[STEPSTRSIZE];
233 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n", alg, ftol, gmx_step_str(count, buf));
235 else if (count < nsteps)
238 "\n%s converged to machine precision in %s steps,\n"
239 "but did not reach the requested Fmax < %g.\n",
241 gmx_step_str(count, buf),
246 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n", alg, ftol, gmx_step_str(count, buf));
250 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
251 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
252 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm / sqrtNumAtoms);
254 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
255 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
256 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm / sqrtNumAtoms);
260 //! Compute the norm and max of the force array in parallel
261 static void get_f_norm_max(const t_commrec* cr,
262 const t_grpopts* opts,
264 gmx::ArrayRef<const gmx::RVec> f,
271 int la_max, a_max, start, end, i, m, gf;
273 /* This routine finds the largest force and returns it.
274 * On parallel machines the global max is taken.
280 end = mdatoms->homenr;
281 if (mdatoms->cFREEZE)
283 for (i = start; i < end; i++)
285 gf = mdatoms->cFREEZE[i];
287 for (m = 0; m < DIM; m++)
289 if (!opts->nFreeze[gf][m])
291 fam += gmx::square(f[i][m]);
304 for (i = start; i < end; i++)
316 if (la_max >= 0 && DOMAINDECOMP(cr))
318 a_max = cr->dd->globalAtomIndices[la_max];
326 snew(sum, 2 * cr->nnodes + 1);
327 sum[2 * cr->nodeid] = fmax2;
328 sum[2 * cr->nodeid + 1] = a_max;
329 sum[2 * cr->nnodes] = fnorm2;
330 gmx_sumd(2 * cr->nnodes + 1, sum, cr);
331 fnorm2 = sum[2 * cr->nnodes];
332 /* Determine the global maximum */
333 for (i = 0; i < cr->nnodes; i++)
335 if (sum[2 * i] > fmax2)
338 a_max = gmx::roundToInt(sum[2 * i + 1]);
346 *fnorm = sqrt(fnorm2);
358 //! Compute the norm of the force
359 static void get_state_f_norm_max(const t_commrec* cr, const t_grpopts* opts, t_mdatoms* mdatoms, em_state_t* ems)
361 get_f_norm_max(cr, opts, mdatoms, ems->f.view().force(), &ems->fnorm, &ems->fmax, &ems->a_fmax);
364 //! Initialize the energy minimization
365 static void init_em(FILE* fplog,
366 const gmx::MDLogger& mdlog,
369 const t_inputrec* ir,
370 gmx::ImdSession* imdSession,
372 t_state* state_global,
373 const gmx_mtop_t& top_global,
378 gmx::MDAtoms* mdAtoms,
379 gmx_global_stat_t* gstat,
380 VirtualSitesHandler* vsite,
381 gmx::Constraints* constr,
382 gmx_shellfc_t** shellfc)
388 fprintf(fplog, "Initiating %s\n", title);
393 state_global->ngtc = 0;
395 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
396 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
397 initialize_lambdas(fplog,
401 ir->simtempvals->temperatures,
402 gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
407 if (ir->eI == IntegrationAlgorithm::NM)
409 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
411 *shellfc = init_shell_flexcon(stdout,
413 constr ? constr->numFlexibleConstraints() : 0,
416 thisRankHasDuty(cr, DUTY_PME));
420 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI),
421 "This else currently only handles energy minimizers, consider if your algorithm "
422 "needs shell/flexible-constraint support");
424 /* With energy minimization, shells and flexible constraints are
425 * automatically minimized when treated like normal DOFS.
427 if (shellfc != nullptr)
433 if (DOMAINDECOMP(cr))
435 dd_init_local_state(*cr->dd, state_global, &ems->s);
437 /* Distribute the charge groups over the nodes from the master node */
438 dd_partition_system(fplog,
459 dd_store_state(*cr->dd, &ems->s);
463 state_change_natoms(state_global, state_global->natoms);
464 /* Just copy the state */
465 ems->s = *state_global;
466 state_change_natoms(&ems->s, ems->s.natoms);
468 mdAlgorithmsSetupAtomData(
469 cr, *ir, top_global, top, fr, &ems->f, mdAtoms, constr, vsite, shellfc ? *shellfc : nullptr);
472 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[FreeEnergyPerturbationCouplingType::Mass]);
476 // TODO how should this cross-module support dependency be managed?
477 if (ir->eConstrAlg == ConstraintAlgorithm::Shake && gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
480 "Can not do energy minimization with %s, use %s\n",
481 enumValueToString(ConstraintAlgorithm::Shake),
482 enumValueToString(ConstraintAlgorithm::Lincs));
485 if (!ir->bContinuation)
487 /* Constrain the starting coordinates */
488 bool needsLogging = true;
489 bool computeEnergy = true;
490 bool computeVirial = false;
492 constr->apply(needsLogging,
497 ems->s.x.arrayRefWithPadding(),
498 ems->s.x.arrayRefWithPadding(),
501 ems->s.lambda[FreeEnergyPerturbationCouplingType::Fep],
503 gmx::ArrayRefWithPadding<RVec>(),
506 gmx::ConstraintVariable::Positions);
512 *gstat = global_stat_init(ir);
519 calc_shifts(ems->s.box, fr->shift_vec);
522 //! Finalize the minimization
523 static void finish_em(const t_commrec* cr,
525 gmx_walltime_accounting_t walltime_accounting,
526 gmx_wallcycle* wcycle)
528 if (!thisRankHasDuty(cr, DUTY_PME))
530 /* Tell the PME only node to finish */
531 gmx_pme_send_finish(cr);
536 em_time_end(walltime_accounting, wcycle);
539 //! Swap two different EM states during minimization
540 static void swap_em_state(em_state_t** ems1, em_state_t** ems2)
549 //! Save the EM trajectory
550 static void write_em_traj(FILE* fplog,
556 const gmx_mtop_t& top_global,
557 const t_inputrec* ir,
560 t_state* state_global,
561 ObservablesHistory* observablesHistory)
567 mdof_flags |= MDOF_X;
571 mdof_flags |= MDOF_F;
574 /* If we want IMD output, set appropriate MDOF flag */
577 mdof_flags |= MDOF_IMD;
580 gmx::WriteCheckpointDataHolder checkpointDataHolder;
581 mdoutf_write_to_trajectory_files(fplog,
587 static_cast<double>(step),
591 state->f.view().force(),
592 &checkpointDataHolder);
594 if (confout != nullptr)
596 if (DOMAINDECOMP(cr))
598 /* If bX=true, x was collected to state_global in the call above */
601 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
603 cr->dd, state->s.ddp_count, state->s.ddp_count_cg_gl, state->s.cg_gl, state->s.x, globalXRef);
608 /* Copy the local state pointer */
609 state_global = &state->s;
614 if (ir->pbcType != PbcType::No && !ir->bPeriodicMols && DOMAINDECOMP(cr))
616 /* Make molecules whole only for confout writing */
617 do_pbc_mtop(ir->pbcType, state->s.box, &top_global, state_global->x.rvec_array());
620 write_sto_conf_mtop(confout,
623 state_global->x.rvec_array(),
631 //! \brief Do one minimization step
633 // \returns true when the step succeeded, false when a constraint error occurred
634 static bool do_em_step(const t_commrec* cr,
635 const t_inputrec* ir,
639 gmx::ArrayRefWithPadding<const gmx::RVec> force,
641 gmx::Constraints* constr,
648 int nthreads gmx_unused;
650 bool validStep = true;
655 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
657 gmx_incons("state mismatch in do_em_step");
660 s2->flags = s1->flags;
662 if (s2->natoms != s1->natoms)
664 state_change_natoms(s2, s1->natoms);
665 ems2->f.resize(s2->natoms);
667 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
669 s2->cg_gl.resize(s1->cg_gl.size());
672 copy_mat(s1->box, s2->box);
673 /* Copy free energy state */
674 s2->lambda = s1->lambda;
675 copy_mat(s1->box, s2->box);
680 nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Update);
681 #pragma omp parallel num_threads(nthreads)
683 const rvec* x1 = s1->x.rvec_array();
684 rvec* x2 = s2->x.rvec_array();
685 const rvec* f = as_rvec_array(force.unpaddedArrayRef().data());
688 #pragma omp for schedule(static) nowait
689 for (int i = start; i < end; i++)
697 for (int m = 0; m < DIM; m++)
699 if (ir->opts.nFreeze[gf][m])
705 x2[i][m] = x1[i][m] + a * f[i][m];
709 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
712 if (s2->flags & enumValueToBitMask(StateEntry::Cgp))
714 /* Copy the CG p vector */
715 const rvec* p1 = s1->cg_p.rvec_array();
716 rvec* p2 = s2->cg_p.rvec_array();
717 #pragma omp for schedule(static) nowait
718 for (int i = start; i < end; i++)
720 // Trivial OpenMP block that does not throw
721 copy_rvec(p1[i], p2[i]);
725 if (DOMAINDECOMP(cr))
727 /* OpenMP does not supported unsigned loop variables */
728 #pragma omp for schedule(static) nowait
729 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
731 s2->cg_gl[i] = s1->cg_gl[i];
736 if (DOMAINDECOMP(cr))
738 s2->ddp_count = s1->ddp_count;
739 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
745 validStep = constr->apply(TRUE,
750 s1->x.arrayRefWithPadding(),
751 s2->x.arrayRefWithPadding(),
754 s2->lambda[FreeEnergyPerturbationCouplingType::Bonded],
756 gmx::ArrayRefWithPadding<RVec>(),
759 gmx::ConstraintVariable::Positions);
763 /* This global reduction will affect performance at high
764 * parallelization, but we can not really avoid it.
765 * But usually EM is not run at high parallelization.
767 int reductionBuffer = static_cast<int>(!validStep);
768 gmx_sumi(1, &reductionBuffer, cr);
769 validStep = (reductionBuffer == 0);
772 // We should move this check to the different minimizers
773 if (!validStep && ir->eI != IntegrationAlgorithm::Steep)
776 "The coordinates could not be constrained. Minimizer '%s' can not handle "
777 "constraint failures, use minimizer '%s' before using '%s'.",
778 enumValueToString(ir->eI),
779 enumValueToString(IntegrationAlgorithm::Steep),
780 enumValueToString(ir->eI));
787 //! Prepare EM for using domain decomposition parallellization
788 static void em_dd_partition_system(FILE* fplog,
789 const gmx::MDLogger& mdlog,
792 const gmx_mtop_t& top_global,
793 const t_inputrec* ir,
794 gmx::ImdSession* imdSession,
798 gmx::MDAtoms* mdAtoms,
800 VirtualSitesHandler* vsite,
801 gmx::Constraints* constr,
803 gmx_wallcycle* wcycle)
805 /* Repartition the domain decomposition */
806 dd_partition_system(fplog,
827 dd_store_state(*cr->dd, &ems->s);
833 //! Copy coordinates, OpenMP parallelized, from \p refCoords to coords
834 void setCoordinates(std::vector<RVec>* coords, ArrayRef<const RVec> refCoords)
836 coords->resize(refCoords.size());
838 const int gmx_unused nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Update);
839 #pragma omp parallel for num_threads(nthreads) schedule(static)
840 for (int i = 0; i < ssize(refCoords); i++)
842 (*coords)[i] = refCoords[i];
846 //! Returns the maximum difference an atom moved between two coordinate sets, over all ranks
847 real maxCoordinateDifference(ArrayRef<const RVec> coords1, ArrayRef<const RVec> coords2, MPI_Comm mpiCommMyGroup)
849 GMX_RELEASE_ASSERT(coords1.size() == coords2.size(), "Coordinate counts should match");
851 real maxDiffSquared = 0;
853 const int gmx_unused nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Update);
854 #pragma omp parallel for reduction(max : maxDiffSquared) num_threads(nthreads) schedule(static)
855 for (int i = 0; i < ssize(coords1); i++)
857 maxDiffSquared = std::max(maxDiffSquared, gmx::norm2(coords1[i] - coords2[i]));
862 if (mpiCommMyGroup != MPI_COMM_NULL)
864 MPI_Comm_size(mpiCommMyGroup, &numRanks);
868 real maxDiffSquaredReduced;
870 &maxDiffSquared, &maxDiffSquaredReduced, 1, GMX_DOUBLE ? MPI_DOUBLE : MPI_FLOAT, MPI_MAX, mpiCommMyGroup);
871 maxDiffSquared = maxDiffSquaredReduced;
874 GMX_UNUSED_VALUE(mpiCommMyGroup);
877 return std::sqrt(maxDiffSquared);
880 /*! \brief Class to handle the work of setting and doing an energy evaluation.
882 * This class is a mere aggregate of parameters to pass to evaluate an
883 * energy, so that future changes to names and types of them consume
884 * less time when refactoring other code.
886 * Aggregate initialization is used, for which the chief risk is that
887 * if a member is added at the end and not all initializer lists are
888 * updated, then the member will be value initialized, which will
889 * typically mean initialization to zero.
891 * Use a braced initializer list to construct one of these. */
892 class EnergyEvaluator
895 /*! \brief Evaluates an energy on the state in \c ems.
897 * \todo In practice, the same objects mu_tot, vir, and pres
898 * are always passed to this function, so we would rather have
899 * them as data members. However, their C-array types are
900 * unsuited for aggregate initialization. When the types
901 * improve, the call signature of this method can be reduced.
903 void run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst);
904 //! Handles logging (deprecated).
907 const gmx::MDLogger& mdlog;
908 //! Handles communication.
910 //! Coordinates multi-simulations.
911 const gmx_multisim_t* ms;
912 //! Holds the simulation topology.
913 const gmx_mtop_t& top_global;
914 //! Holds the domain topology.
916 //! User input options.
917 const t_inputrec* inputrec;
918 //! The Interactive Molecular Dynamics session.
919 gmx::ImdSession* imdSession;
920 //! The pull work object.
922 //! Manages flop accounting.
924 //! Manages wall cycle accounting.
925 gmx_wallcycle* wcycle;
926 //! Coordinates global reduction.
927 gmx_global_stat_t gstat;
928 //! Handles virtual sites.
929 VirtualSitesHandler* vsite;
930 //! Handles constraints.
931 gmx::Constraints* constr;
932 //! Per-atom data for this domain.
933 gmx::MDAtoms* mdAtoms;
934 //! Handles how to calculate the forces.
936 //! Schedule of force-calculation work each step for this task.
937 MdrunScheduleWorkload* runScheduleWork;
938 //! Stores the computed energies.
939 gmx_enerdata_t* enerd;
940 //! The DD partitioning count at which the pair list was generated
941 int ddpCountPairSearch;
942 //! The local coordinates that were used for pair searching, stored for computing displacements
943 std::vector<RVec> pairSearchCoordinates;
946 void EnergyEvaluator::run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst)
950 tensor force_vir, shake_vir, ekin;
954 /* Set the time to the initial time, the time does not change during EM */
955 t = inputrec->init_t;
959 vsite->construct(ems->s.x, {}, ems->s.box, gmx::VSiteOperation::Positions);
962 // Compute the buffer size of the pair list
963 const real bufferSize = inputrec->rlist - std::max(inputrec->rcoulomb, inputrec->rvdw);
965 if (bFirst || bufferSize <= 0 || (DOMAINDECOMP(cr) && ems->s.ddp_count != ddpCountPairSearch))
967 /* This is the first state or an old state used before the last ns */
972 // We need to generate a new pairlist when one atom moved more than half the buffer size
973 ArrayRef<const RVec> localCoordinates =
974 ArrayRef<const RVec>(ems->s.x).subArray(0, mdAtoms->mdatoms()->homenr);
975 bNS = 2 * maxCoordinateDifference(pairSearchCoordinates, localCoordinates, cr->mpi_comm_mygroup)
979 if (DOMAINDECOMP(cr) && bNS)
981 /* Repartition the domain decomposition */
982 em_dd_partition_system(
983 fplog, mdlog, count, cr, top_global, inputrec, imdSession, pull_work, ems, top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
984 ddpCountPairSearch = cr->dd->ddp_count;
987 /* Store the local coordinates that will be used in the pair search, after we re-partitioned */
988 if (bufferSize > 0 && bNS)
990 ArrayRef<const RVec> localCoordinates =
991 constArrayRefFromArray(ems->s.x.data(), mdAtoms->mdatoms()->homenr);
992 setCoordinates(&pairSearchCoordinates, localCoordinates);
995 /* Calc force & energy on new trial position */
996 /* do_force always puts the charge groups in the box and shifts again
997 * We do not unshift, so molecules are always whole in congrad.c
1012 ems->s.x.arrayRefWithPadding(),
1025 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY
1026 | (bNS ? GMX_FORCE_NS : 0),
1027 DDBalanceRegionHandler(cr));
1029 /* Clear the unused shake virial and pressure */
1030 clear_mat(shake_vir);
1033 /* Communicate stuff when parallel */
1034 if (PAR(cr) && inputrec->eI != IntegrationAlgorithm::NM)
1036 wallcycle_start(wcycle, WallCycleCounter::MoveE);
1045 gmx::ArrayRef<real>{},
1047 std::vector<real>(1, terminate),
1049 CGLO_ENERGY | CGLO_PRESSURE | CGLO_CONSTRAINT);
1051 wallcycle_stop(wcycle, WallCycleCounter::MoveE);
1054 if (fr->dispersionCorrection)
1056 /* Calculate long range corrections to pressure and energy */
1057 const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
1058 ems->s.box, ems->s.lambda[FreeEnergyPerturbationCouplingType::Vdw]);
1060 enerd->term[F_DISPCORR] = correction.energy;
1061 enerd->term[F_EPOT] += correction.energy;
1062 enerd->term[F_PRES] += correction.pressure;
1063 enerd->term[F_DVDL] += correction.dvdl;
1067 enerd->term[F_DISPCORR] = 0;
1070 ems->epot = enerd->term[F_EPOT];
1074 /* Project out the constraint components of the force */
1075 bool needsLogging = false;
1076 bool computeEnergy = false;
1077 bool computeVirial = true;
1079 auto f = ems->f.view().forceWithPadding();
1080 constr->apply(needsLogging,
1085 ems->s.x.arrayRefWithPadding(),
1087 f.unpaddedArrayRef(),
1089 ems->s.lambda[FreeEnergyPerturbationCouplingType::Bonded],
1091 gmx::ArrayRefWithPadding<RVec>(),
1094 gmx::ConstraintVariable::ForceDispl);
1095 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1096 m_add(force_vir, shake_vir, vir);
1100 copy_mat(force_vir, vir);
1104 enerd->term[F_PRES] = calc_pres(fr->pbcType, inputrec->nwall, ems->s.box, ekin, vir, pres);
1106 if (inputrec->efep != FreeEnergyPerturbationType::No)
1108 accumulateKineticLambdaComponents(enerd, ems->s.lambda, *inputrec->fepvals);
1111 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
1113 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
1119 //! Parallel utility summing energies and forces
1120 static double reorder_partsum(const t_commrec* cr,
1121 const t_grpopts* opts,
1122 const gmx_mtop_t& top_global,
1123 const em_state_t* s_min,
1124 const em_state_t* s_b)
1128 fprintf(debug, "Doing reorder_partsum\n");
1131 auto fm = s_min->f.view().force();
1132 auto fb = s_b->f.view().force();
1134 /* Collect fm in a global vector fmg.
1135 * This conflicts with the spirit of domain decomposition,
1136 * but to fully optimize this a much more complicated algorithm is required.
1138 const int natoms = top_global.natoms;
1142 gmx::ArrayRef<const int> indicesMin = s_min->s.cg_gl;
1144 for (int a : indicesMin)
1146 copy_rvec(fm[i], fmg[a]);
1149 gmx_sum(top_global.natoms * 3, fmg[0], cr);
1151 /* Now we will determine the part of the sum for the cgs in state s_b */
1152 gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
1157 gmx::ArrayRef<const unsigned char> grpnrFREEZE =
1158 top_global.groups.groupNumbers[SimulationAtomGroupType::Freeze];
1159 for (int a : indicesB)
1161 if (!grpnrFREEZE.empty())
1163 gf = grpnrFREEZE[i];
1165 for (int m = 0; m < DIM; m++)
1167 if (!opts->nFreeze[gf][m])
1169 partsum += (fb[i][m] - fmg[a][m]) * fb[i][m];
1180 //! Print some stuff, like beta, whatever that means.
1181 static real pr_beta(const t_commrec* cr,
1182 const t_grpopts* opts,
1184 const gmx_mtop_t& top_global,
1185 const em_state_t* s_min,
1186 const em_state_t* s_b)
1190 /* This is just the classical Polak-Ribiere calculation of beta;
1191 * it looks a bit complicated since we take freeze groups into account,
1192 * and might have to sum it in parallel runs.
1195 if (!DOMAINDECOMP(cr)
1196 || (s_min->s.ddp_count == cr->dd->ddp_count && s_b->s.ddp_count == cr->dd->ddp_count))
1198 auto fm = s_min->f.view().force();
1199 auto fb = s_b->f.view().force();
1202 /* This part of code can be incorrect with DD,
1203 * since the atom ordering in s_b and s_min might differ.
1205 for (int i = 0; i < mdatoms->homenr; i++)
1207 if (mdatoms->cFREEZE)
1209 gf = mdatoms->cFREEZE[i];
1211 for (int m = 0; m < DIM; m++)
1213 if (!opts->nFreeze[gf][m])
1215 sum += (fb[i][m] - fm[i][m]) * fb[i][m];
1222 /* We need to reorder cgs while summing */
1223 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1227 gmx_sumd(1, &sum, cr);
1230 return sum / gmx::square(s_min->fnorm);
1236 void LegacySimulator::do_cg()
1238 const char* CG = "Polak-Ribiere Conjugate Gradients";
1240 gmx_localtop_t top(top_global.ffparams);
1241 gmx_global_stat_t gstat;
1242 double tmp, minstep;
1244 real a, b, c, beta = 0.0;
1247 gmx_bool converged, foundlower;
1248 rvec mu_tot = { 0 };
1249 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1251 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1252 int m, step, nminstep;
1253 auto* mdatoms = mdAtoms->mdatoms();
1258 "Note that activating conjugate gradient energy minimization via the "
1259 "integrator .mdp option and the command gmx mdrun may "
1260 "be available in a different form in a future version of GROMACS, "
1261 "e.g. gmx minimize and an .mdp option.");
1267 // In CG, the state is extended with a search direction
1268 state_global->flags |= enumValueToBitMask(StateEntry::Cgp);
1270 // Ensure the extra per-atom state array gets allocated
1271 state_change_natoms(state_global, state_global->natoms);
1273 // Initialize the search direction to zero
1274 for (RVec& cg_p : state_global->cg_p)
1280 /* Create 4 states on the stack and extract pointers that we will swap */
1281 em_state_t s0{}, s1{}, s2{}, s3{};
1282 em_state_t* s_min = &s0;
1283 em_state_t* s_a = &s1;
1284 em_state_t* s_b = &s2;
1285 em_state_t* s_c = &s3;
1287 /* Init em and store the local state in s_min */
1306 const bool simulationsShareState = false;
1307 gmx_mdoutf* outf = init_mdoutf(fplog,
1318 StartingBehavior::NewSimulation,
1319 simulationsShareState,
1321 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
1327 StartingBehavior::NewSimulation,
1328 simulationsShareState,
1329 mdModulesNotifiers);
1331 /* Print to log file */
1332 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1334 /* Max number of steps */
1335 number_steps = inputrec->nsteps;
1339 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1343 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1346 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global,
1347 &top, inputrec, imdSession, pull_work, nrnb,
1348 wcycle, gstat, vsite, constr, mdAtoms,
1349 fr, runScheduleWork, enerd, -1, {} };
1350 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1351 /* do_force always puts the charge groups in the box and shifts again
1352 * We do not unshift, so molecules are always whole in congrad.c
1354 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1358 /* Copy stuff to the energy bin for easy printing etc. */
1359 matrix nullBox = {};
1360 energyOutput.addDataAtEnergyStep(false,
1362 static_cast<double>(step),
1378 EnergyOutput::printHeader(fplog, step, step);
1379 energyOutput.printStepToEnergyFile(
1380 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
1383 /* Estimate/guess the initial stepsize */
1384 stepsize = inputrec->em_stepsize / s_min->fnorm;
1388 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1389 fprintf(stderr, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1390 fprintf(stderr, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1391 fprintf(stderr, "\n");
1392 /* and copy to the log file too... */
1393 fprintf(fplog, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1394 fprintf(fplog, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1395 fprintf(fplog, "\n");
1397 /* Start the loop over CG steps.
1398 * Each successful step is counted, and we continue until
1399 * we either converge or reach the max number of steps.
1402 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1405 /* start taking steps in a new direction
1406 * First time we enter the routine, beta=0, and the direction is
1407 * simply the negative gradient.
1410 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1411 gmx::ArrayRef<gmx::RVec> pm = s_min->s.cg_p;
1412 gmx::ArrayRef<const gmx::RVec> sfm = s_min->f.view().force();
1415 for (int i = 0; i < mdatoms->homenr; i++)
1417 if (mdatoms->cFREEZE)
1419 gf = mdatoms->cFREEZE[i];
1421 for (m = 0; m < DIM; m++)
1423 if (!inputrec->opts.nFreeze[gf][m])
1425 pm[i][m] = sfm[i][m] + beta * pm[i][m];
1426 gpa -= pm[i][m] * sfm[i][m];
1427 /* f is negative gradient, thus the sign */
1436 /* Sum the gradient along the line across CPUs */
1439 gmx_sumd(1, &gpa, cr);
1442 /* Calculate the norm of the search vector */
1443 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1445 /* Just in case stepsize reaches zero due to numerical precision... */
1448 stepsize = inputrec->em_stepsize / pnorm;
1452 * Double check the value of the derivative in the search direction.
1453 * If it is positive it must be due to the old information in the
1454 * CG formula, so just remove that and start over with beta=0.
1455 * This corresponds to a steepest descent step.
1460 step--; /* Don't count this step since we are restarting */
1461 continue; /* Go back to the beginning of the big for-loop */
1464 /* Calculate minimum allowed stepsize, before the average (norm)
1465 * relative change in coordinate is smaller than precision
1468 auto s_min_x = makeArrayRef(s_min->s.x);
1469 for (int i = 0; i < mdatoms->homenr; i++)
1471 for (m = 0; m < DIM; m++)
1473 tmp = fabs(s_min_x[i][m]);
1478 tmp = pm[i][m] / tmp;
1479 minstep += tmp * tmp;
1482 /* Add up from all CPUs */
1485 gmx_sumd(1, &minstep, cr);
1488 minstep = GMX_REAL_EPS / sqrt(minstep / (3 * top_global.natoms));
1490 if (stepsize < minstep)
1496 /* Write coordinates if necessary */
1497 do_x = do_per_step(step, inputrec->nstxout);
1498 do_f = do_per_step(step, inputrec->nstfout);
1501 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, step, s_min, state_global, observablesHistory);
1503 /* Take a step downhill.
1504 * In theory, we should minimize the function along this direction.
1505 * That is quite possible, but it turns out to take 5-10 function evaluations
1506 * for each line. However, we dont really need to find the exact minimum -
1507 * it is much better to start a new CG step in a modified direction as soon
1508 * as we are close to it. This will save a lot of energy evaluations.
1510 * In practice, we just try to take a single step.
1511 * If it worked (i.e. lowered the energy), we increase the stepsize but
1512 * the continue straight to the next CG step without trying to find any minimum.
1513 * If it didn't work (higher energy), there must be a minimum somewhere between
1514 * the old position and the new one.
1516 * Due to the finite numerical accuracy, it turns out that it is a good idea
1517 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1518 * This leads to lower final energies in the tests I've done. / Erik
1520 s_a->epot = s_min->epot;
1522 c = a + stepsize; /* reference position along line is zero */
1524 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1526 em_dd_partition_system(fplog,
1544 /* Take a trial step (new coords in s_c) */
1545 do_em_step(cr, inputrec, mdatoms, s_min, c, s_min->s.cg_p.constArrayRefWithPadding(), s_c, constr, -1);
1548 /* Calculate energy for the trial step */
1549 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1551 /* Calc derivative along line */
1552 const rvec* pc = s_c->s.cg_p.rvec_array();
1553 gmx::ArrayRef<const gmx::RVec> sfc = s_c->f.view().force();
1555 for (int i = 0; i < mdatoms->homenr; i++)
1557 for (m = 0; m < DIM; m++)
1559 gpc -= pc[i][m] * sfc[i][m]; /* f is negative gradient, thus the sign */
1562 /* Sum the gradient along the line across CPUs */
1565 gmx_sumd(1, &gpc, cr);
1568 /* This is the max amount of increase in energy we tolerate */
1569 tmp = std::sqrt(GMX_REAL_EPS) * fabs(s_a->epot);
1571 /* Accept the step if the energy is lower, or if it is not significantly higher
1572 * and the line derivative is still negative.
1574 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1577 /* Great, we found a better energy. Increase step for next iteration
1578 * if we are still going down, decrease it otherwise
1582 stepsize *= 1.618034; /* The golden section */
1586 stepsize *= 0.618034; /* 1/golden section */
1591 /* New energy is the same or higher. We will have to do some work
1592 * to find a smaller value in the interval. Take smaller step next time!
1595 stepsize *= 0.618034;
1599 /* OK, if we didn't find a lower value we will have to locate one now - there must
1600 * be one in the interval [a=0,c].
1601 * The same thing is valid here, though: Don't spend dozens of iterations to find
1602 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1603 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1605 * I also have a safeguard for potentially really pathological functions so we never
1606 * take more than 20 steps before we give up ...
1608 * If we already found a lower value we just skip this step and continue to the update.
1617 /* Select a new trial point.
1618 * If the derivatives at points a & c have different sign we interpolate to zero,
1619 * otherwise just do a bisection.
1621 if (gpa < 0 && gpc > 0)
1623 b = a + gpa * (a - c) / (gpc - gpa);
1630 /* safeguard if interpolation close to machine accuracy causes errors:
1631 * never go outside the interval
1633 if (b <= a || b >= c)
1638 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1640 /* Reload the old state */
1641 em_dd_partition_system(fplog,
1659 /* Take a trial step to this new point - new coords in s_b */
1660 do_em_step(cr, inputrec, mdatoms, s_min, b, s_min->s.cg_p.constArrayRefWithPadding(), s_b, constr, -1);
1663 /* Calculate energy for the trial step */
1664 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1666 /* p does not change within a step, but since the domain decomposition
1667 * might change, we have to use cg_p of s_b here.
1669 const rvec* pb = s_b->s.cg_p.rvec_array();
1670 gmx::ArrayRef<const gmx::RVec> sfb = s_b->f.view().force();
1672 for (int i = 0; i < mdatoms->homenr; i++)
1674 for (m = 0; m < DIM; m++)
1676 gpb -= pb[i][m] * sfb[i][m]; /* f is negative gradient, thus the sign */
1679 /* Sum the gradient along the line across CPUs */
1682 gmx_sumd(1, &gpb, cr);
1687 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a->epot, s_b->epot, s_c->epot, gpb);
1690 epot_repl = s_b->epot;
1692 /* Keep one of the intervals based on the value of the derivative at the new point */
1695 /* Replace c endpoint with b */
1696 swap_em_state(&s_b, &s_c);
1702 /* Replace a endpoint with b */
1703 swap_em_state(&s_b, &s_a);
1709 * Stop search as soon as we find a value smaller than the endpoints.
1710 * Never run more than 20 steps, no matter what.
1713 } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) && (nminstep < 20));
1715 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot) * GMX_REAL_EPS || nminstep >= 20)
1717 /* OK. We couldn't find a significantly lower energy.
1718 * If beta==0 this was steepest descent, and then we give up.
1719 * If not, set beta=0 and restart with steepest descent before quitting.
1729 /* Reset memory before giving up */
1735 /* Select min energy state of A & C, put the best in B.
1737 if (s_c->epot < s_a->epot)
1741 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c->epot, s_a->epot);
1743 swap_em_state(&s_b, &s_c);
1750 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a->epot, s_c->epot);
1752 swap_em_state(&s_b, &s_a);
1760 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n", s_c->epot);
1762 swap_em_state(&s_b, &s_c);
1766 /* new search direction */
1767 /* beta = 0 means forget all memory and restart with steepest descents. */
1768 if (nstcg && ((step % nstcg) == 0))
1774 /* s_min->fnorm cannot be zero, because then we would have converged
1778 /* Polak-Ribiere update.
1779 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1781 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1783 /* Limit beta to prevent oscillations */
1784 if (fabs(beta) > 5.0)
1790 /* update positions */
1791 swap_em_state(&s_min, &s_b);
1794 /* Print it if necessary */
1797 if (mdrunOptions.verbose)
1799 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1801 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1804 s_min->fnorm / sqrtNumAtoms,
1809 /* Store the new (lower) energies */
1810 matrix nullBox = {};
1811 energyOutput.addDataAtEnergyStep(false,
1813 static_cast<double>(step),
1829 do_log = do_per_step(step, inputrec->nstlog);
1830 do_ene = do_per_step(step, inputrec->nstenergy);
1832 imdSession->fillEnergyRecord(step, TRUE);
1836 EnergyOutput::printHeader(fplog, step, step);
1838 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1842 do_log ? fplog : nullptr,
1849 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1850 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x, 0))
1852 imdSession->sendPositionsAndEnergies();
1855 /* Stop when the maximum force lies below tolerance.
1856 * If we have reached machine precision, converged is already set to true.
1858 converged = converged || (s_min->fmax < inputrec->em_tol);
1860 } /* End of the loop */
1864 step--; /* we never took that last step in this case */
1866 if (s_min->fmax > inputrec->em_tol)
1870 warn_step(fplog, inputrec->em_tol, s_min->fmax, step - 1 == number_steps, FALSE);
1877 /* If we printed energy and/or logfile last step (which was the last step)
1878 * we don't have to do it again, but otherwise print the final values.
1882 /* Write final value to log since we didn't do anything the last step */
1883 EnergyOutput::printHeader(fplog, step, step);
1885 if (!do_ene || !do_log)
1887 /* Write final energy file entries */
1888 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1892 !do_log ? fplog : nullptr,
1900 /* Print some stuff... */
1903 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1907 * For accurate normal mode calculation it is imperative that we
1908 * store the last conformation into the full precision binary trajectory.
1910 * However, we should only do it if we did NOT already write this step
1911 * above (which we did if do_x or do_f was true).
1913 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1914 * in the trajectory with the same step number.
1916 do_x = !do_per_step(step, inputrec->nstxout);
1917 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1920 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, s_min, state_global, observablesHistory);
1925 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1926 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1927 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1929 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1932 finish_em(cr, outf, walltime_accounting, wcycle);
1934 /* To print the actual number of steps we needed somewhere */
1935 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1939 void LegacySimulator::do_lbfgs()
1941 static const char* LBFGS = "Low-Memory BFGS Minimizer";
1943 gmx_localtop_t top(top_global.ffparams);
1944 gmx_global_stat_t gstat;
1945 auto* mdatoms = mdAtoms->mdatoms();
1950 "Note that activating L-BFGS energy minimization via the "
1951 "integrator .mdp option and the command gmx mdrun may "
1952 "be available in a different form in a future version of GROMACS, "
1953 "e.g. gmx minimize and an .mdp option.");
1957 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1960 if (nullptr != constr)
1964 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1965 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1968 const int n = 3 * state_global->natoms;
1969 const int nmaxcorr = inputrec->nbfgscorr;
1971 std::vector<real> p(n);
1972 std::vector<real> rho(nmaxcorr);
1973 std::vector<real> alpha(nmaxcorr);
1975 std::vector<std::vector<real>> dx(nmaxcorr);
1976 for (auto& dxCorr : dx)
1981 std::vector<std::vector<real>> dg(nmaxcorr);
1982 for (auto& dgCorr : dg)
2009 const bool simulationsShareState = false;
2010 gmx_mdoutf* outf = init_mdoutf(fplog,
2021 StartingBehavior::NewSimulation,
2022 simulationsShareState,
2024 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
2030 StartingBehavior::NewSimulation,
2031 simulationsShareState,
2032 mdModulesNotifiers);
2034 const int start = 0;
2035 const int end = mdatoms->homenr;
2037 /* We need 4 working states */
2038 em_state_t s0{}, s1{}, s2{}, s3{};
2039 em_state_t* sa = &s0;
2040 em_state_t* sb = &s1;
2041 em_state_t* sc = &s2;
2042 em_state_t* last = &s3;
2043 /* Initialize by copying the state from ems (we could skip x and f here) */
2048 /* Print to log file */
2049 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
2051 /* Max number of steps */
2052 const int number_steps = inputrec->nsteps;
2054 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
2055 std::vector<bool> frozen(n);
2057 for (int i = start; i < end; i++)
2059 if (mdatoms->cFREEZE)
2061 gf = mdatoms->cFREEZE[i];
2063 for (int m = 0; m < DIM; m++)
2065 frozen[3 * i + m] = (inputrec->opts.nFreeze[gf][m] != 0);
2070 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
2074 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
2079 vsite->construct(state_global->x, {}, state_global->box, VSiteOperation::Positions);
2082 /* Call the force routine and some auxiliary (neighboursearching etc.) */
2083 /* do_force always puts the charge groups in the box and shifts again
2084 * We do not unshift, so molecules are always whole
2087 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2088 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2089 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2093 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
2097 /* Copy stuff to the energy bin for easy printing etc. */
2098 matrix nullBox = {};
2099 energyOutput.addDataAtEnergyStep(false,
2101 static_cast<double>(step),
2117 EnergyOutput::printHeader(fplog, step, step);
2118 energyOutput.printStepToEnergyFile(
2119 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
2122 /* Set the initial step.
2123 * since it will be multiplied by the non-normalized search direction
2124 * vector (force vector the first time), we scale it by the
2125 * norm of the force.
2130 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2131 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2132 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2133 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2134 fprintf(stderr, "\n");
2135 /* and copy to the log file too... */
2136 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2137 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2138 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2139 fprintf(fplog, "\n");
2142 // Point is an index to the memory of search directions, where 0 is the first one.
2145 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
2146 real* fInit = static_cast<real*>(ems.f.view().force().data()[0]);
2147 for (int i = 0; i < n; i++)
2151 dx[point][i] = fInit[i]; /* Initial search direction */
2159 // Stepsize will be modified during the search, and actually it is not critical
2160 // (the main efficiency in the algorithm comes from changing directions), but
2161 // we still need an initial value, so estimate it as the inverse of the norm
2162 // so we take small steps where the potential fluctuates a lot.
2163 double stepsize = 1.0 / ems.fnorm;
2165 /* Start the loop over BFGS steps.
2166 * Each successful step is counted, and we continue until
2167 * we either converge or reach the max number of steps.
2175 /* Set the gradient from the force */
2176 bool converged = false;
2177 for (int step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
2180 /* Write coordinates if necessary */
2181 const bool do_x = do_per_step(step, inputrec->nstxout);
2182 const bool do_f = do_per_step(step, inputrec->nstfout);
2187 mdof_flags |= MDOF_X;
2192 mdof_flags |= MDOF_F;
2197 mdof_flags |= MDOF_IMD;
2200 gmx::WriteCheckpointDataHolder checkpointDataHolder;
2201 mdoutf_write_to_trajectory_files(fplog,
2207 static_cast<real>(step),
2211 ems.f.view().force(),
2212 &checkpointDataHolder);
2214 /* Do the linesearching in the direction dx[point][0..(n-1)] */
2216 /* make s a pointer to current search direction - point=0 first time we get here */
2217 gmx::ArrayRef<const real> s = dx[point];
2219 const real* xx = static_cast<real*>(ems.s.x.rvec_array()[0]);
2220 const real* ff = static_cast<real*>(ems.f.view().force().data()[0]);
2222 // calculate line gradient in position A
2224 for (int i = 0; i < n; i++)
2226 gpa -= s[i] * ff[i];
2229 /* Calculate minimum allowed stepsize along the line, before the average (norm)
2230 * relative change in coordinate is smaller than precision
2233 for (int i = 0; i < n; i++)
2235 double tmp = fabs(xx[i]);
2241 minstep += tmp * tmp;
2243 minstep = GMX_REAL_EPS / sqrt(minstep / n);
2245 if (stepsize < minstep)
2251 // Before taking any steps along the line, store the old position
2253 real* lastx = static_cast<real*>(last->s.x.data()[0]);
2254 real* lastf = static_cast<real*>(last->f.view().force().data()[0]);
2255 const real Epot0 = ems.epot;
2259 /* Take a step downhill.
2260 * In theory, we should find the actual minimum of the function in this
2261 * direction, somewhere along the line.
2262 * That is quite possible, but it turns out to take 5-10 function evaluations
2263 * for each line. However, we dont really need to find the exact minimum -
2264 * it is much better to start a new BFGS step in a modified direction as soon
2265 * as we are close to it. This will save a lot of energy evaluations.
2267 * In practice, we just try to take a single step.
2268 * If it worked (i.e. lowered the energy), we increase the stepsize but
2269 * continue straight to the next BFGS step without trying to find any minimum,
2270 * i.e. we change the search direction too. If the line was smooth, it is
2271 * likely we are in a smooth region, and then it makes sense to take longer
2272 * steps in the modified search direction too.
2274 * If it didn't work (higher energy), there must be a minimum somewhere between
2275 * the old position and the new one. Then we need to start by finding a lower
2276 * value before we change search direction. Since the energy was apparently
2277 * quite rough, we need to decrease the step size.
2279 * Due to the finite numerical accuracy, it turns out that it is a good idea
2280 * to accept a SMALL increase in energy, if the derivative is still downhill.
2281 * This leads to lower final energies in the tests I've done. / Erik
2284 // State "A" is the first position along the line.
2285 // reference position along line is initially zero
2288 // Check stepsize first. We do not allow displacements
2289 // larger than emstep.
2295 // Pick a new position C by adding stepsize to A.
2298 // Calculate what the largest change in any individual coordinate
2299 // would be (translation along line * gradient along line)
2301 for (int i = 0; i < n; i++)
2303 real delta = c * s[i];
2304 if (delta > maxdelta)
2309 // If any displacement is larger than the stepsize limit, reduce the step
2310 if (maxdelta > inputrec->em_stepsize)
2314 } while (maxdelta > inputrec->em_stepsize);
2316 // Take a trial step and move the coordinate array xc[] to position C
2317 real* xc = static_cast<real*>(sc->s.x.rvec_array()[0]);
2318 for (int i = 0; i < n; i++)
2320 xc[i] = lastx[i] + c * s[i];
2324 // Calculate energy for the trial step in position C
2325 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2327 // Calc line gradient in position C
2328 real* fc = static_cast<real*>(sc->f.view().force()[0]);
2330 for (int i = 0; i < n; i++)
2332 gpc -= s[i] * fc[i]; /* f is negative gradient, thus the sign */
2334 /* Sum the gradient along the line across CPUs */
2337 gmx_sumd(1, &gpc, cr);
2340 // This is the max amount of increase in energy we tolerate.
2341 // By allowing VERY small changes (close to numerical precision) we
2342 // frequently find even better (lower) final energies.
2343 double tmp = std::sqrt(GMX_REAL_EPS) * fabs(sa->epot);
2345 // Accept the step if the energy is lower in the new position C (compared to A),
2346 // or if it is not significantly higher and the line derivative is still negative.
2347 bool foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2348 // If true, great, we found a better energy. We no longer try to alter the
2349 // stepsize, but simply accept this new better position. The we select a new
2350 // search direction instead, which will be much more efficient than continuing
2351 // to take smaller steps along a line. Set fnorm based on the new C position,
2352 // which will be used to update the stepsize to 1/fnorm further down.
2354 // If false, the energy is NOT lower in point C, i.e. it will be the same
2355 // or higher than in point A. In this case it is pointless to move to point C,
2356 // so we will have to do more iterations along the same line to find a smaller
2357 // value in the interval [A=0.0,C].
2358 // Here, A is still 0.0, but that will change when we do a search in the interval
2359 // [0.0,C] below. That search we will do by interpolation or bisection rather
2360 // than with the stepsize, so no need to modify it. For the next search direction
2361 // it will be reset to 1/fnorm anyway.
2366 // OK, if we didn't find a lower value we will have to locate one now - there must
2367 // be one in the interval [a,c].
2368 // The same thing is valid here, though: Don't spend dozens of iterations to find
2369 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2370 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2371 // I also have a safeguard for potentially really pathological functions so we never
2372 // take more than 20 steps before we give up.
2373 // If we already found a lower value we just skip this step and continue to the update.
2378 // Select a new trial point B in the interval [A,C].
2379 // If the derivatives at points a & c have different sign we interpolate to zero,
2380 // otherwise just do a bisection since there might be multiple minima/maxima
2381 // inside the interval.
2383 if (gpa < 0 && gpc > 0)
2385 b = a + gpa * (a - c) / (gpc - gpa);
2392 /* safeguard if interpolation close to machine accuracy causes errors:
2393 * never go outside the interval
2395 if (b <= a || b >= c)
2400 // Take a trial step to point B
2401 real* xb = static_cast<real*>(sb->s.x.rvec_array()[0]);
2402 for (int i = 0; i < n; i++)
2404 xb[i] = lastx[i] + b * s[i];
2408 // Calculate energy for the trial step in point B
2409 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2412 // Calculate gradient in point B
2413 real* fb = static_cast<real*>(sb->f.view().force()[0]);
2415 for (int i = 0; i < n; i++)
2417 gpb -= s[i] * fb[i]; /* f is negative gradient, thus the sign */
2419 /* Sum the gradient along the line across CPUs */
2422 gmx_sumd(1, &gpb, cr);
2425 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2426 // at the new point B, and rename the endpoints of this new interval A and C.
2429 /* Replace c endpoint with b */
2431 /* copy state b to c */
2436 /* Replace a endpoint with b */
2438 /* copy state b to a */
2443 * Stop search as soon as we find a value smaller than the endpoints,
2444 * or if the tolerance is below machine precision.
2445 * Never run more than 20 steps, no matter what.
2448 } while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2450 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2452 /* OK. We couldn't find a significantly lower energy.
2453 * If ncorr==0 this was steepest descent, and then we give up.
2454 * If not, reset memory to restart as steepest descent before quitting.
2466 /* Search in gradient direction */
2467 for (int i = 0; i < n; i++)
2469 dx[point][i] = ff[i];
2471 /* Reset stepsize */
2472 stepsize = 1.0 / fnorm;
2477 /* Select min energy state of A & C, put the best in xx/ff/Epot
2479 if (sc->epot < sa->epot)
2500 /* Update the memory information, and calculate a new
2501 * approximation of the inverse hessian
2504 /* Have new data in Epot, xx, ff */
2505 if (ncorr < nmaxcorr)
2510 for (int i = 0; i < n; i++)
2512 dg[point][i] = lastf[i] - ff[i];
2513 dx[point][i] *= step_taken;
2518 for (int i = 0; i < n; i++)
2520 dgdg += dg[point][i] * dg[point][i];
2521 dgdx += dg[point][i] * dx[point][i];
2524 const real diag = dgdx / dgdg;
2526 rho[point] = 1.0 / dgdx;
2529 if (point >= nmaxcorr)
2535 for (int i = 0; i < n; i++)
2542 /* Recursive update. First go back over the memory points */
2543 for (int k = 0; k < ncorr; k++)
2552 for (int i = 0; i < n; i++)
2554 sq += dx[cp][i] * p[i];
2557 alpha[cp] = rho[cp] * sq;
2559 for (int i = 0; i < n; i++)
2561 p[i] -= alpha[cp] * dg[cp][i];
2565 for (int i = 0; i < n; i++)
2570 /* And then go forward again */
2571 for (int k = 0; k < ncorr; k++)
2574 for (int i = 0; i < n; i++)
2576 yr += p[i] * dg[cp][i];
2579 real beta = rho[cp] * yr;
2580 beta = alpha[cp] - beta;
2582 for (int i = 0; i < n; i++)
2584 p[i] += beta * dx[cp][i];
2594 for (int i = 0; i < n; i++)
2598 dx[point][i] = p[i];
2606 /* Print it if necessary */
2609 if (mdrunOptions.verbose)
2611 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2613 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2616 ems.fnorm / sqrtNumAtoms,
2621 /* Store the new (lower) energies */
2622 matrix nullBox = {};
2623 energyOutput.addDataAtEnergyStep(false,
2625 static_cast<double>(step),
2641 do_log = do_per_step(step, inputrec->nstlog);
2642 do_ene = do_per_step(step, inputrec->nstenergy);
2644 imdSession->fillEnergyRecord(step, TRUE);
2648 EnergyOutput::printHeader(fplog, step, step);
2650 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2654 do_log ? fplog : nullptr,
2661 /* Send x and E to IMD client, if bIMD is TRUE. */
2662 if (imdSession->run(step, TRUE, state_global->box, state_global->x, 0) && MASTER(cr))
2664 imdSession->sendPositionsAndEnergies();
2667 // Reset stepsize in we are doing more iterations
2670 /* Stop when the maximum force lies below tolerance.
2671 * If we have reached machine precision, converged is already set to true.
2673 converged = converged || (ems.fmax < inputrec->em_tol);
2675 } /* End of the loop */
2679 step--; /* we never took that last step in this case */
2681 if (ems.fmax > inputrec->em_tol)
2685 warn_step(fplog, inputrec->em_tol, ems.fmax, step - 1 == number_steps, FALSE);
2690 /* If we printed energy and/or logfile last step (which was the last step)
2691 * we don't have to do it again, but otherwise print the final values.
2693 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2695 EnergyOutput::printHeader(fplog, step, step);
2697 if (!do_ene || !do_log) /* Write final energy file entries */
2699 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2703 !do_log ? fplog : nullptr,
2710 /* Print some stuff... */
2713 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2717 * For accurate normal mode calculation it is imperative that we
2718 * store the last conformation into the full precision binary trajectory.
2720 * However, we should only do it if we did NOT already write this step
2721 * above (which we did if do_x or do_f was true).
2723 const bool do_x = !do_per_step(step, inputrec->nstxout);
2724 const bool do_f = !do_per_step(step, inputrec->nstfout);
2726 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, &ems, state_global, observablesHistory);
2730 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2731 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2732 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2734 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2737 finish_em(cr, outf, walltime_accounting, wcycle);
2739 /* To print the actual number of steps we needed somewhere */
2740 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2743 void LegacySimulator::do_steep()
2745 const char* SD = "Steepest Descents";
2746 gmx_localtop_t top(top_global.ffparams);
2747 gmx_global_stat_t gstat;
2750 gmx_bool bDone, bAbort, do_x, do_f;
2752 rvec mu_tot = { 0 };
2755 int steps_accepted = 0;
2756 auto* mdatoms = mdAtoms->mdatoms();
2761 "Note that activating steepest-descent energy minimization via the "
2762 "integrator .mdp option and the command gmx mdrun may "
2763 "be available in a different form in a future version of GROMACS, "
2764 "e.g. gmx minimize and an .mdp option.");
2766 /* Create 2 states on the stack and extract pointers that we will swap */
2767 em_state_t s0{}, s1{};
2768 em_state_t* s_min = &s0;
2769 em_state_t* s_try = &s1;
2771 /* Init em and store the local state in s_try */
2790 const bool simulationsShareState = false;
2791 gmx_mdoutf* outf = init_mdoutf(fplog,
2802 StartingBehavior::NewSimulation,
2803 simulationsShareState,
2805 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
2811 StartingBehavior::NewSimulation,
2812 simulationsShareState,
2813 mdModulesNotifiers);
2815 /* Print to log file */
2816 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2818 /* Set variables for stepsize (in nm). This is the largest
2819 * step that we are going to make in any direction.
2821 ustep = inputrec->em_stepsize;
2824 /* Max number of steps */
2825 nsteps = inputrec->nsteps;
2829 /* Print to the screen */
2830 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2834 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2836 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2837 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2838 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2840 /**** HERE STARTS THE LOOP ****
2841 * count is the counter for the number of steps
2842 * bDone will be TRUE when the minimization has converged
2843 * bAbort will be TRUE when nsteps steps have been performed or when
2844 * the stepsize becomes smaller than is reasonable for machine precision
2849 while (!bDone && !bAbort)
2851 bAbort = (nsteps >= 0) && (count == nsteps);
2853 /* set new coordinates, except for first step */
2854 bool validStep = true;
2857 validStep = do_em_step(
2858 cr, inputrec, mdatoms, s_min, stepsize, s_min->f.view().forceWithPadding(), s_try, constr, count);
2863 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2867 // Signal constraint error during stepping with energy=inf
2868 s_try->epot = std::numeric_limits<real>::infinity();
2873 EnergyOutput::printHeader(fplog, count, count);
2878 s_min->epot = s_try->epot;
2881 /* Print it if necessary */
2884 if (mdrunOptions.verbose)
2887 "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2893 ((count == 0) || (s_try->epot < s_min->epot)) ? '\n' : '\r');
2897 if ((count == 0) || (s_try->epot < s_min->epot))
2899 /* Store the new (lower) energies */
2900 matrix nullBox = {};
2901 energyOutput.addDataAtEnergyStep(false,
2903 static_cast<double>(count),
2919 imdSession->fillEnergyRecord(count, TRUE);
2921 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2922 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2923 energyOutput.printStepToEnergyFile(
2924 mdoutf_get_fp_ene(outf), TRUE, do_dr, do_or, fplog, count, count, fr->fcdata.get(), nullptr);
2929 /* Now if the new energy is smaller than the previous...
2930 * or if this is the first step!
2931 * or if we did random steps!
2934 if ((count == 0) || (s_try->epot < s_min->epot))
2938 /* Test whether the convergence criterion is met... */
2939 bDone = (s_try->fmax < inputrec->em_tol);
2941 /* Copy the arrays for force, positions and energy */
2942 /* The 'Min' array always holds the coords and forces of the minimal
2944 swap_em_state(&s_min, &s_try);
2950 /* Write to trn, if necessary */
2951 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2952 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2954 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, count, s_min, state_global, observablesHistory);
2958 /* If energy is not smaller make the step smaller... */
2961 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2963 /* Reload the old state */
2964 em_dd_partition_system(fplog,
2983 // If the force is very small after finishing minimization,
2984 // we risk dividing by zero when calculating the step size.
2985 // So we check first if the minimization has stopped before
2986 // trying to obtain a new step size.
2989 /* Determine new step */
2990 stepsize = ustep / s_min->fmax;
2993 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2995 if (count == nsteps || ustep < 1e-12)
2997 if (count == nsteps || ustep < 1e-6)
3002 warn_step(fplog, inputrec->em_tol, s_min->fmax, count == nsteps, constr != nullptr);
3007 /* Send IMD energies and positions, if bIMD is TRUE. */
3008 if (imdSession->run(count,
3010 MASTER(cr) ? state_global->box : nullptr,
3011 MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>(),
3015 imdSession->sendPositionsAndEnergies();
3019 } /* End of the loop */
3021 /* Print some data... */
3024 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
3026 write_em_traj(fplog,
3030 inputrec->nstfout != 0,
3031 ftp2fn(efSTO, nfile, fnm),
3037 observablesHistory);
3041 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
3043 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
3044 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
3047 finish_em(cr, outf, walltime_accounting, wcycle);
3049 /* To print the actual number of steps we needed somewhere */
3051 // TODO: Avoid changing inputrec (#3854)
3052 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
3053 nonConstInputrec->nsteps = count;
3056 walltime_accounting_set_nsteps_done(walltime_accounting, count);
3059 void LegacySimulator::do_nm()
3061 const char* NM = "Normal Mode Analysis";
3063 gmx_localtop_t top(top_global.ffparams);
3064 gmx_global_stat_t gstat;
3066 rvec mu_tot = { 0 };
3068 gmx_bool bSparse; /* use sparse matrix storage format */
3070 gmx_sparsematrix_t* sparse_matrix = nullptr;
3071 real* full_matrix = nullptr;
3073 /* added with respect to mdrun */
3075 real der_range = 10.0 * std::sqrt(GMX_REAL_EPS);
3077 bool bIsMaster = MASTER(cr);
3078 auto* mdatoms = mdAtoms->mdatoms();
3083 "Note that activating normal-mode analysis via the integrator "
3084 ".mdp option and the command gmx mdrun may "
3085 "be available in a different form in a future version of GROMACS, "
3086 "e.g. gmx normal-modes.");
3088 if (constr != nullptr)
3092 "Constraints present with Normal Mode Analysis, this combination is not supported");
3095 gmx_shellfc_t* shellfc;
3097 em_state_t state_work{};
3099 /* Init em and store the local state in state_minimum */
3118 const bool simulationsShareState = false;
3119 gmx_mdoutf* outf = init_mdoutf(fplog,
3130 StartingBehavior::NewSimulation,
3131 simulationsShareState,
3134 std::vector<int> atom_index = get_atom_index(top_global);
3135 std::vector<gmx::RVec> fneg(atom_index.size(), { 0, 0, 0 });
3136 snew(dfdx, atom_index.size());
3142 "NOTE: This version of GROMACS has been compiled in single precision,\n"
3143 " which MIGHT not be accurate enough for normal mode analysis.\n"
3144 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
3145 " are fairly modest even if you recompile in double precision.\n\n");
3149 /* Check if we can/should use sparse storage format.
3151 * Sparse format is only useful when the Hessian itself is sparse, which it
3152 * will be when we use a cutoff.
3153 * For small systems (n<1000) it is easier to always use full matrix format, though.
3155 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
3157 GMX_LOG(mdlog.warning)
3158 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
3161 else if (atom_index.size() < 1000)
3163 GMX_LOG(mdlog.warning)
3164 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
3170 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
3174 /* Number of dimensions, based on real atoms, that is not vsites or shell */
3175 sz = DIM * atom_index.size();
3177 fprintf(stderr, "Allocating Hessian memory...\n\n");
3181 sparse_matrix = gmx_sparsematrix_init(sz);
3182 sparse_matrix->compressed_symmetric = TRUE;
3186 snew(full_matrix, sz * sz);
3189 /* Write start time and temperature */
3190 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
3192 /* fudge nr of steps to nr of atoms */
3194 // TODO: Avoid changing inputrec (#3854)
3195 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
3196 nonConstInputrec->nsteps = atom_index.size() * 2;
3202 "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
3207 nnodes = cr->nnodes;
3209 /* Make evaluate_energy do a single node force calculation */
3211 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
3212 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
3213 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
3214 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
3215 cr->nnodes = nnodes;
3217 /* if forces are not small, warn user */
3218 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
3220 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
3221 if (state_work.fmax > 1.0e-3)
3223 GMX_LOG(mdlog.warning)
3225 "The force is probably not small enough to "
3226 "ensure that you are at a minimum.\n"
3227 "Be aware that negative eigenvalues may occur\n"
3228 "when the resulting matrix is diagonalized.");
3231 /***********************************************************
3233 * Loop over all pairs in matrix
3235 * do_force called twice. Once with positive and
3236 * once with negative displacement
3238 ************************************************************/
3240 /* Steps are divided one by one over the nodes */
3242 auto state_work_x = makeArrayRef(state_work.s.x);
3243 auto state_work_f = state_work.f.view().force();
3244 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
3246 size_t atom = atom_index[aid];
3247 for (size_t d = 0; d < DIM; d++)
3250 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
3253 x_min = state_work_x[atom][d];
3255 for (unsigned int dx = 0; (dx < 2); dx++)
3259 state_work_x[atom][d] = x_min - der_range;
3263 state_work_x[atom][d] = x_min + der_range;
3266 /* Make evaluate_energy do a single node force calculation */
3270 /* Now is the time to relax the shells */
3271 relax_shell_flexcon(fplog,
3274 mdrunOptions.verbose,
3285 state_work.s.natoms,
3286 state_work.s.x.arrayRefWithPadding(),
3287 state_work.s.v.arrayRefWithPadding(),
3289 state_work.s.lambda,
3291 &state_work.f.view(),
3302 DDBalanceRegionHandler(nullptr));
3308 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid * 2 + dx, FALSE);
3311 cr->nnodes = nnodes;
3315 std::copy(state_work_f.begin(), state_work_f.begin() + atom_index.size(), fneg.begin());
3319 /* x is restored to original */
3320 state_work_x[atom][d] = x_min;
3322 for (size_t j = 0; j < atom_index.size(); j++)
3324 for (size_t k = 0; (k < DIM); k++)
3326 dfdx[j][k] = -(state_work_f[atom_index[j]][k] - fneg[j][k]) / (2 * der_range);
3333 # define mpi_type GMX_MPI_REAL
3334 MPI_Send(dfdx[0], atom_index.size() * DIM, mpi_type, MASTER(cr), cr->nodeid, cr->mpi_comm_mygroup);
3339 for (index node = 0; (node < nnodes && aid + node < ssize(atom_index)); node++)
3345 MPI_Recv(dfdx[0], atom_index.size() * DIM, mpi_type, node, node, cr->mpi_comm_mygroup, &stat);
3350 row = (aid + node) * DIM + d;
3352 for (size_t j = 0; j < atom_index.size(); j++)
3354 for (size_t k = 0; k < DIM; k++)
3360 if (col >= row && dfdx[j][k] != 0.0)
3362 gmx_sparsematrix_increment_value(sparse_matrix, row, col, dfdx[j][k]);
3367 full_matrix[row * sz + col] = dfdx[j][k];
3374 if (mdrunOptions.verbose && fplog)
3379 /* write progress */
3380 if (bIsMaster && mdrunOptions.verbose)
3383 "\rFinished step %d out of %td",
3384 std::min<int>(atom + nnodes, atom_index.size()),
3392 fprintf(stderr, "\n\nWriting Hessian...\n");
3393 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
3396 finish_em(cr, outf, walltime_accounting, wcycle);
3398 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size() * 2);