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
58 #include "gromacs/commandline/filenm.h"
59 #include "gromacs/domdec/collect.h"
60 #include "gromacs/domdec/dlbtiming.h"
61 #include "gromacs/domdec/domdec.h"
62 #include "gromacs/domdec/domdec_struct.h"
63 #include "gromacs/domdec/mdsetup.h"
64 #include "gromacs/domdec/partition.h"
65 #include "gromacs/ewald/pme_pp.h"
66 #include "gromacs/fileio/confio.h"
67 #include "gromacs/fileio/mtxio.h"
68 #include "gromacs/gmxlib/network.h"
69 #include "gromacs/gmxlib/nrnb.h"
70 #include "gromacs/imd/imd.h"
71 #include "gromacs/linearalgebra/sparsematrix.h"
72 #include "gromacs/listed_forces/listed_forces.h"
73 #include "gromacs/math/functions.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/mdlib/constr.h"
76 #include "gromacs/mdlib/coupling.h"
77 #include "gromacs/mdlib/dispersioncorrection.h"
78 #include "gromacs/mdlib/ebin.h"
79 #include "gromacs/mdlib/enerdata_utils.h"
80 #include "gromacs/mdlib/energyoutput.h"
81 #include "gromacs/mdlib/force.h"
82 #include "gromacs/mdlib/force_flags.h"
83 #include "gromacs/mdlib/forcerec.h"
84 #include "gromacs/mdlib/gmx_omp_nthreads.h"
85 #include "gromacs/mdlib/md_support.h"
86 #include "gromacs/mdlib/mdatoms.h"
87 #include "gromacs/mdlib/stat.h"
88 #include "gromacs/mdlib/tgroup.h"
89 #include "gromacs/mdlib/trajectory_writing.h"
90 #include "gromacs/mdlib/update.h"
91 #include "gromacs/mdlib/vsite.h"
92 #include "gromacs/mdrunutility/handlerestart.h"
93 #include "gromacs/mdrunutility/printtime.h"
94 #include "gromacs/mdtypes/checkpointdata.h"
95 #include "gromacs/mdtypes/commrec.h"
96 #include "gromacs/mdtypes/forcebuffers.h"
97 #include "gromacs/mdtypes/forcerec.h"
98 #include "gromacs/mdtypes/inputrec.h"
99 #include "gromacs/mdtypes/interaction_const.h"
100 #include "gromacs/mdtypes/md_enums.h"
101 #include "gromacs/mdtypes/mdatom.h"
102 #include "gromacs/mdtypes/mdrunoptions.h"
103 #include "gromacs/mdtypes/observablesreducer.h"
104 #include "gromacs/mdtypes/state.h"
105 #include "gromacs/pbcutil/pbc.h"
106 #include "gromacs/timing/wallcycle.h"
107 #include "gromacs/timing/walltime_accounting.h"
108 #include "gromacs/topology/mtop_util.h"
109 #include "gromacs/topology/topology.h"
110 #include "gromacs/utility/cstringutil.h"
111 #include "gromacs/utility/exceptions.h"
112 #include "gromacs/utility/fatalerror.h"
113 #include "gromacs/utility/logger.h"
114 #include "gromacs/utility/smalloc.h"
116 #include "legacysimulator.h"
120 using gmx::MdrunScheduleWorkload;
122 using gmx::VirtualSitesHandler;
124 //! Utility structure for manipulating states during EM
125 typedef struct em_state
127 //! Copy of the global state
133 //! Norm of the force
141 //! Print the EM starting conditions
142 static void print_em_start(FILE* fplog,
144 gmx_walltime_accounting_t walltime_accounting,
145 gmx_wallcycle* wcycle,
148 walltime_accounting_start_time(walltime_accounting);
149 wallcycle_start(wcycle, WallCycleCounter::Run);
150 print_start(fplog, cr, walltime_accounting, name);
153 //! Stop counting time for EM
154 static void em_time_end(gmx_walltime_accounting_t walltime_accounting, gmx_wallcycle* wcycle)
156 wallcycle_stop(wcycle, WallCycleCounter::Run);
158 walltime_accounting_end_time(walltime_accounting);
161 //! Printing a log file and console header
162 static void sp_header(FILE* out, const char* minimizer, real ftol, int nsteps)
165 fprintf(out, "%s:\n", minimizer);
166 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
167 fprintf(out, " Number of steps = %12d\n", nsteps);
170 //! Print warning message
171 static void warn_step(FILE* fp, real ftol, real fmax, gmx_bool bLastStep, gmx_bool bConstrain)
173 constexpr bool realIsDouble = GMX_DOUBLE;
176 if (!std::isfinite(fmax))
179 "\nEnergy minimization has stopped because the force "
180 "on at least one atom is not finite. This usually means "
181 "atoms are overlapping. Modify the input coordinates to "
182 "remove atom overlap or use soft-core potentials with "
183 "the free energy code to avoid infinite forces.\n%s",
184 !realIsDouble ? "You could also be lucky that switching to double precision "
185 "is sufficient to obtain finite forces.\n"
191 "\nEnergy minimization reached the maximum number "
192 "of steps before the forces reached the requested "
193 "precision Fmax < %g.\n",
199 "\nEnergy minimization has stopped, but the forces have "
200 "not converged to the requested precision Fmax < %g (which "
201 "may not be possible for your system). It stopped "
202 "because the algorithm tried to make a new step whose size "
203 "was too small, or there was no change in the energy since "
204 "last step. Either way, we regard the minimization as "
205 "converged to within the available machine precision, "
206 "given your starting configuration and EM parameters.\n%s%s",
208 !realIsDouble ? "\nDouble precision normally gives you higher accuracy, but "
209 "this is often not needed for preparing to run molecular "
212 bConstrain ? "You might need to increase your constraint accuracy, or turn\n"
213 "off constraints altogether (set constraints = none in mdp file)\n"
217 fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
218 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
221 //! Print message about convergence of the EM
222 static void print_converged(FILE* fp,
228 const em_state_t* ems,
231 char buf[STEPSTRSIZE];
235 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n", alg, ftol, gmx_step_str(count, buf));
237 else if (count < nsteps)
240 "\n%s converged to machine precision in %s steps,\n"
241 "but did not reach the requested Fmax < %g.\n",
243 gmx_step_str(count, buf),
248 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n", alg, ftol, gmx_step_str(count, buf));
252 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
253 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
254 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm / sqrtNumAtoms);
256 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
257 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
258 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm / sqrtNumAtoms);
262 //! Compute the norm and max of the force array in parallel
263 static void get_f_norm_max(const t_commrec* cr,
264 const t_grpopts* opts,
266 gmx::ArrayRef<const gmx::RVec> f,
273 int la_max, a_max, start, end, i, m, gf;
275 /* This routine finds the largest force and returns it.
276 * On parallel machines the global max is taken.
282 end = mdatoms->homenr;
283 if (mdatoms->cFREEZE)
285 for (i = start; i < end; i++)
287 gf = mdatoms->cFREEZE[i];
289 for (m = 0; m < DIM; m++)
291 if (!opts->nFreeze[gf][m])
293 fam += gmx::square(f[i][m]);
306 for (i = start; i < end; i++)
318 if (la_max >= 0 && DOMAINDECOMP(cr))
320 a_max = cr->dd->globalAtomIndices[la_max];
328 snew(sum, 2 * cr->nnodes + 1);
329 sum[2 * cr->nodeid] = fmax2;
330 sum[2 * cr->nodeid + 1] = a_max;
331 sum[2 * cr->nnodes] = fnorm2;
332 gmx_sumd(2 * cr->nnodes + 1, sum, cr);
333 fnorm2 = sum[2 * cr->nnodes];
334 /* Determine the global maximum */
335 for (i = 0; i < cr->nnodes; i++)
337 if (sum[2 * i] > fmax2)
340 a_max = gmx::roundToInt(sum[2 * i + 1]);
348 *fnorm = sqrt(fnorm2);
360 //! Compute the norm of the force
361 static void get_state_f_norm_max(const t_commrec* cr, const t_grpopts* opts, t_mdatoms* mdatoms, em_state_t* ems)
363 get_f_norm_max(cr, opts, mdatoms, ems->f.view().force(), &ems->fnorm, &ems->fmax, &ems->a_fmax);
366 //! Initialize the energy minimization
367 static void init_em(FILE* fplog,
368 const gmx::MDLogger& mdlog,
371 const t_inputrec* ir,
372 gmx::ImdSession* imdSession,
374 t_state* state_global,
375 const gmx_mtop_t& top_global,
380 gmx::MDAtoms* mdAtoms,
381 gmx_global_stat_t* gstat,
382 VirtualSitesHandler* vsite,
383 gmx::Constraints* constr,
384 gmx_shellfc_t** shellfc)
390 fprintf(fplog, "Initiating %s\n", title);
395 state_global->ngtc = 0;
397 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
398 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
399 initialize_lambdas(fplog,
403 ir->simtempvals->temperatures,
404 gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc),
409 if (ir->eI == IntegrationAlgorithm::NM)
411 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
413 *shellfc = init_shell_flexcon(stdout,
415 constr ? constr->numFlexibleConstraints() : 0,
418 thisRankHasDuty(cr, DUTY_PME));
422 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI),
423 "This else currently only handles energy minimizers, consider if your algorithm "
424 "needs shell/flexible-constraint support");
426 /* With energy minimization, shells and flexible constraints are
427 * automatically minimized when treated like normal DOFS.
429 if (shellfc != nullptr)
435 if (DOMAINDECOMP(cr))
437 // Local state only becomes valid now.
438 dd_init_local_state(*cr->dd, state_global, &ems->s);
440 /* Distribute the charge groups over the nodes from the master node */
441 dd_partition_system(fplog,
462 dd_store_state(*cr->dd, &ems->s);
466 state_change_natoms(state_global, state_global->natoms);
467 /* Just copy the state */
468 ems->s = *state_global;
469 state_change_natoms(&ems->s, ems->s.natoms);
471 mdAlgorithmsSetupAtomData(
472 cr, *ir, top_global, top, fr, &ems->f, mdAtoms, constr, vsite, shellfc ? *shellfc : nullptr);
475 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[FreeEnergyPerturbationCouplingType::Mass]);
479 // TODO how should this cross-module support dependency be managed?
480 if (ir->eConstrAlg == ConstraintAlgorithm::Shake && gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
483 "Can not do energy minimization with %s, use %s\n",
484 enumValueToString(ConstraintAlgorithm::Shake),
485 enumValueToString(ConstraintAlgorithm::Lincs));
488 if (!ir->bContinuation)
490 /* Constrain the starting coordinates */
491 bool needsLogging = true;
492 bool computeEnergy = true;
493 bool computeVirial = false;
495 constr->apply(needsLogging,
500 ems->s.x.arrayRefWithPadding(),
501 ems->s.x.arrayRefWithPadding(),
504 ems->s.lambda[FreeEnergyPerturbationCouplingType::Fep],
506 gmx::ArrayRefWithPadding<RVec>(),
509 gmx::ConstraintVariable::Positions);
515 *gstat = global_stat_init(ir);
522 calc_shifts(ems->s.box, fr->shift_vec);
525 //! Finalize the minimization
526 static void finish_em(const t_commrec* cr,
528 gmx_walltime_accounting_t walltime_accounting,
529 gmx_wallcycle* wcycle)
531 if (!thisRankHasDuty(cr, DUTY_PME))
533 /* Tell the PME only node to finish */
534 gmx_pme_send_finish(cr);
539 em_time_end(walltime_accounting, wcycle);
542 //! Swap two different EM states during minimization
543 static void swap_em_state(em_state_t** ems1, em_state_t** ems2)
552 //! Save the EM trajectory
553 static void write_em_traj(FILE* fplog,
559 const gmx_mtop_t& top_global,
560 const t_inputrec* ir,
563 t_state* state_global,
564 ObservablesHistory* observablesHistory)
570 mdof_flags |= MDOF_X;
574 mdof_flags |= MDOF_F;
577 /* If we want IMD output, set appropriate MDOF flag */
580 mdof_flags |= MDOF_IMD;
583 gmx::WriteCheckpointDataHolder checkpointDataHolder;
584 mdoutf_write_to_trajectory_files(fplog,
590 static_cast<double>(step),
594 state->f.view().force(),
595 &checkpointDataHolder);
597 if (confout != nullptr)
599 if (DOMAINDECOMP(cr))
601 /* If bX=true, x was collected to state_global in the call above */
604 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
606 cr->dd, state->s.ddp_count, state->s.ddp_count_cg_gl, state->s.cg_gl, state->s.x, globalXRef);
611 /* Copy the local state pointer */
612 state_global = &state->s;
617 if (ir->pbcType != PbcType::No && !ir->bPeriodicMols && DOMAINDECOMP(cr))
619 /* Make molecules whole only for confout writing */
620 do_pbc_mtop(ir->pbcType, state->s.box, &top_global, state_global->x.rvec_array());
623 write_sto_conf_mtop(confout,
626 state_global->x.rvec_array(),
634 //! \brief Do one minimization step
636 // \returns true when the step succeeded, false when a constraint error occurred
637 static bool do_em_step(const t_commrec* cr,
638 const t_inputrec* ir,
642 gmx::ArrayRefWithPadding<const gmx::RVec> force,
644 gmx::Constraints* constr,
651 int nthreads gmx_unused;
653 bool validStep = true;
658 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
660 gmx_incons("state mismatch in do_em_step");
663 s2->flags = s1->flags;
665 if (s2->natoms != s1->natoms)
667 state_change_natoms(s2, s1->natoms);
668 ems2->f.resize(s2->natoms);
670 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
672 s2->cg_gl.resize(s1->cg_gl.size());
675 copy_mat(s1->box, s2->box);
676 /* Copy free energy state */
677 s2->lambda = s1->lambda;
678 copy_mat(s1->box, s2->box);
683 nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Update);
684 #pragma omp parallel num_threads(nthreads)
686 const rvec* x1 = s1->x.rvec_array();
687 rvec* x2 = s2->x.rvec_array();
688 const rvec* f = as_rvec_array(force.unpaddedArrayRef().data());
691 #pragma omp for schedule(static) nowait
692 for (int i = start; i < end; i++)
700 for (int m = 0; m < DIM; m++)
702 if (ir->opts.nFreeze[gf][m])
708 x2[i][m] = x1[i][m] + a * f[i][m];
712 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
715 if (s2->flags & enumValueToBitMask(StateEntry::Cgp))
717 /* Copy the CG p vector */
718 const rvec* p1 = s1->cg_p.rvec_array();
719 rvec* p2 = s2->cg_p.rvec_array();
720 #pragma omp for schedule(static) nowait
721 for (int i = start; i < end; i++)
723 // Trivial OpenMP block that does not throw
724 copy_rvec(p1[i], p2[i]);
728 if (DOMAINDECOMP(cr))
730 /* OpenMP does not supported unsigned loop variables */
731 #pragma omp for schedule(static) nowait
732 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
734 s2->cg_gl[i] = s1->cg_gl[i];
739 if (DOMAINDECOMP(cr))
741 s2->ddp_count = s1->ddp_count;
742 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
748 validStep = constr->apply(TRUE,
753 s1->x.arrayRefWithPadding(),
754 s2->x.arrayRefWithPadding(),
757 s2->lambda[FreeEnergyPerturbationCouplingType::Bonded],
759 gmx::ArrayRefWithPadding<RVec>(),
762 gmx::ConstraintVariable::Positions);
766 /* This global reduction will affect performance at high
767 * parallelization, but we can not really avoid it.
768 * But usually EM is not run at high parallelization.
770 int reductionBuffer = static_cast<int>(!validStep);
771 gmx_sumi(1, &reductionBuffer, cr);
772 validStep = (reductionBuffer == 0);
775 // We should move this check to the different minimizers
776 if (!validStep && ir->eI != IntegrationAlgorithm::Steep)
779 "The coordinates could not be constrained. Minimizer '%s' can not handle "
780 "constraint failures, use minimizer '%s' before using '%s'.",
781 enumValueToString(ir->eI),
782 enumValueToString(IntegrationAlgorithm::Steep),
783 enumValueToString(ir->eI));
790 //! Prepare EM for using domain decomposition parallellization
791 static void em_dd_partition_system(FILE* fplog,
792 const gmx::MDLogger& mdlog,
795 const gmx_mtop_t& top_global,
796 const t_inputrec* ir,
797 gmx::ImdSession* imdSession,
801 gmx::MDAtoms* mdAtoms,
803 VirtualSitesHandler* vsite,
804 gmx::Constraints* constr,
806 gmx_wallcycle* wcycle)
808 /* Repartition the domain decomposition */
809 dd_partition_system(fplog,
830 dd_store_state(*cr->dd, &ems->s);
836 //! Copy coordinates, OpenMP parallelized, from \p refCoords to coords
837 void setCoordinates(std::vector<RVec>* coords, ArrayRef<const RVec> refCoords)
839 coords->resize(refCoords.size());
841 const int gmx_unused nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Update);
842 #pragma omp parallel for num_threads(nthreads) schedule(static)
843 for (int i = 0; i < ssize(refCoords); i++)
845 (*coords)[i] = refCoords[i];
849 //! Returns the maximum difference an atom moved between two coordinate sets, over all ranks
850 real maxCoordinateDifference(ArrayRef<const RVec> coords1, ArrayRef<const RVec> coords2, MPI_Comm mpiCommMyGroup)
852 GMX_RELEASE_ASSERT(coords1.size() == coords2.size(), "Coordinate counts should match");
854 real maxDiffSquared = 0;
856 const int gmx_unused nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Update);
857 #pragma omp parallel for reduction(max : maxDiffSquared) num_threads(nthreads) schedule(static)
858 for (int i = 0; i < ssize(coords1); i++)
860 maxDiffSquared = std::max(maxDiffSquared, gmx::norm2(coords1[i] - coords2[i]));
865 if (mpiCommMyGroup != MPI_COMM_NULL)
867 MPI_Comm_size(mpiCommMyGroup, &numRanks);
871 real maxDiffSquaredReduced;
873 &maxDiffSquared, &maxDiffSquaredReduced, 1, GMX_DOUBLE ? MPI_DOUBLE : MPI_FLOAT, MPI_MAX, mpiCommMyGroup);
874 maxDiffSquared = maxDiffSquaredReduced;
877 GMX_UNUSED_VALUE(mpiCommMyGroup);
880 return std::sqrt(maxDiffSquared);
883 /*! \brief Class to handle the work of setting and doing an energy evaluation.
885 * This class is a mere aggregate of parameters to pass to evaluate an
886 * energy, so that future changes to names and types of them consume
887 * less time when refactoring other code.
889 * Aggregate initialization is used, for which the chief risk is that
890 * if a member is added at the end and not all initializer lists are
891 * updated, then the member will be value initialized, which will
892 * typically mean initialization to zero.
894 * Use a braced initializer list to construct one of these. */
895 class EnergyEvaluator
898 /*! \brief Evaluates an energy on the state in \c ems.
900 * \todo In practice, the same objects mu_tot, vir, and pres
901 * are always passed to this function, so we would rather have
902 * them as data members. However, their C-array types are
903 * unsuited for aggregate initialization. When the types
904 * improve, the call signature of this method can be reduced.
906 void run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst, int64_t step);
907 //! Handles logging (deprecated).
910 const gmx::MDLogger& mdlog;
911 //! Handles communication.
913 //! Coordinates multi-simulations.
914 const gmx_multisim_t* ms;
915 //! Holds the simulation topology.
916 const gmx_mtop_t& top_global;
917 //! Holds the domain topology.
919 //! User input options.
920 const t_inputrec* inputrec;
921 //! The Interactive Molecular Dynamics session.
922 gmx::ImdSession* imdSession;
923 //! The pull work object.
925 //! Manages flop accounting.
927 //! Manages wall cycle accounting.
928 gmx_wallcycle* wcycle;
929 //! Legacy coordinator of global reduction.
930 gmx_global_stat_t gstat;
931 //! Coordinates reduction for observables
932 gmx::ObservablesReducer* observablesReducer;
933 //! Handles virtual sites.
934 VirtualSitesHandler* vsite;
935 //! Handles constraints.
936 gmx::Constraints* constr;
937 //! Per-atom data for this domain.
938 gmx::MDAtoms* mdAtoms;
939 //! Handles how to calculate the forces.
941 //! Schedule of force-calculation work each step for this task.
942 MdrunScheduleWorkload* runScheduleWork;
943 //! Stores the computed energies.
944 gmx_enerdata_t* enerd;
945 //! The DD partitioning count at which the pair list was generated
946 int ddpCountPairSearch;
947 //! The local coordinates that were used for pair searching, stored for computing displacements
948 std::vector<RVec> pairSearchCoordinates;
951 void EnergyEvaluator::run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst, int64_t step)
955 tensor force_vir, shake_vir, ekin;
959 /* Set the time to the initial time, the time does not change during EM */
960 t = inputrec->init_t;
964 vsite->construct(ems->s.x, {}, ems->s.box, gmx::VSiteOperation::Positions);
967 // Compute the buffer size of the pair list
968 const real bufferSize = inputrec->rlist - std::max(inputrec->rcoulomb, inputrec->rvdw);
970 if (bFirst || bufferSize <= 0 || (DOMAINDECOMP(cr) && ems->s.ddp_count != ddpCountPairSearch))
972 /* This is the first state or an old state used before the last ns */
977 // We need to generate a new pairlist when one atom moved more than half the buffer size
978 ArrayRef<const RVec> localCoordinates =
979 ArrayRef<const RVec>(ems->s.x).subArray(0, mdAtoms->mdatoms()->homenr);
980 bNS = 2 * maxCoordinateDifference(pairSearchCoordinates, localCoordinates, cr->mpi_comm_mygroup)
984 if (DOMAINDECOMP(cr) && bNS)
986 /* Repartition the domain decomposition */
987 em_dd_partition_system(
988 fplog, mdlog, count, cr, top_global, inputrec, imdSession, pull_work, ems, top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
989 ddpCountPairSearch = cr->dd->ddp_count;
992 /* Store the local coordinates that will be used in the pair search, after we re-partitioned */
993 if (bufferSize > 0 && bNS)
995 ArrayRef<const RVec> localCoordinates =
996 constArrayRefFromArray(ems->s.x.data(), mdAtoms->mdatoms()->homenr);
997 setCoordinates(&pairSearchCoordinates, localCoordinates);
1000 fr->longRangeNonbondeds->updateAfterPartition(*mdAtoms->mdatoms());
1002 /* Calc force & energy on new trial position */
1003 /* do_force always puts the charge groups in the box and shifts again
1004 * We do not unshift, so molecules are always whole in congrad.c
1019 ems->s.x.arrayRefWithPadding(),
1032 fr->longRangeNonbondeds.get(),
1033 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY
1034 | (bNS ? GMX_FORCE_NS : 0),
1035 DDBalanceRegionHandler(cr));
1037 /* Clear the unused shake virial and pressure */
1038 clear_mat(shake_vir);
1041 /* Communicate stuff when parallel */
1042 if (PAR(cr) && inputrec->eI != IntegrationAlgorithm::NM)
1044 wallcycle_start(wcycle, WallCycleCounter::MoveE);
1053 gmx::ArrayRef<real>{},
1055 std::vector<real>(1, terminate),
1057 CGLO_ENERGY | CGLO_PRESSURE | CGLO_CONSTRAINT,
1059 observablesReducer);
1061 wallcycle_stop(wcycle, WallCycleCounter::MoveE);
1064 if (fr->dispersionCorrection)
1066 /* Calculate long range corrections to pressure and energy */
1067 const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
1068 ems->s.box, ems->s.lambda[FreeEnergyPerturbationCouplingType::Vdw]);
1070 enerd->term[F_DISPCORR] = correction.energy;
1071 enerd->term[F_EPOT] += correction.energy;
1072 enerd->term[F_PRES] += correction.pressure;
1073 enerd->term[F_DVDL] += correction.dvdl;
1077 enerd->term[F_DISPCORR] = 0;
1080 ems->epot = enerd->term[F_EPOT];
1084 /* Project out the constraint components of the force */
1085 bool needsLogging = false;
1086 bool computeEnergy = false;
1087 bool computeVirial = true;
1089 auto f = ems->f.view().forceWithPadding();
1090 constr->apply(needsLogging,
1095 ems->s.x.arrayRefWithPadding(),
1097 f.unpaddedArrayRef(),
1099 ems->s.lambda[FreeEnergyPerturbationCouplingType::Bonded],
1101 gmx::ArrayRefWithPadding<RVec>(),
1104 gmx::ConstraintVariable::ForceDispl);
1105 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1106 m_add(force_vir, shake_vir, vir);
1110 copy_mat(force_vir, vir);
1114 enerd->term[F_PRES] = calc_pres(fr->pbcType, inputrec->nwall, ems->s.box, ekin, vir, pres);
1116 if (inputrec->efep != FreeEnergyPerturbationType::No)
1118 accumulateKineticLambdaComponents(enerd, ems->s.lambda, *inputrec->fepvals);
1121 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
1123 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
1129 //! Parallel utility summing energies and forces
1130 static double reorder_partsum(const t_commrec* cr,
1131 const t_grpopts* opts,
1132 const gmx_mtop_t& top_global,
1133 const em_state_t* s_min,
1134 const em_state_t* s_b)
1138 fprintf(debug, "Doing reorder_partsum\n");
1141 auto fm = s_min->f.view().force();
1142 auto fb = s_b->f.view().force();
1144 /* Collect fm in a global vector fmg.
1145 * This conflicts with the spirit of domain decomposition,
1146 * but to fully optimize this a much more complicated algorithm is required.
1148 const int natoms = top_global.natoms;
1152 gmx::ArrayRef<const int> indicesMin = s_min->s.cg_gl;
1154 for (int a : indicesMin)
1156 copy_rvec(fm[i], fmg[a]);
1159 gmx_sum(top_global.natoms * 3, fmg[0], cr);
1161 /* Now we will determine the part of the sum for the cgs in state s_b */
1162 gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
1167 gmx::ArrayRef<const unsigned char> grpnrFREEZE =
1168 top_global.groups.groupNumbers[SimulationAtomGroupType::Freeze];
1169 for (int a : indicesB)
1171 if (!grpnrFREEZE.empty())
1173 gf = grpnrFREEZE[i];
1175 for (int m = 0; m < DIM; m++)
1177 if (!opts->nFreeze[gf][m])
1179 partsum += (fb[i][m] - fmg[a][m]) * fb[i][m];
1190 //! Print some stuff, like beta, whatever that means.
1191 static real pr_beta(const t_commrec* cr,
1192 const t_grpopts* opts,
1194 const gmx_mtop_t& top_global,
1195 const em_state_t* s_min,
1196 const em_state_t* s_b)
1200 /* This is just the classical Polak-Ribiere calculation of beta;
1201 * it looks a bit complicated since we take freeze groups into account,
1202 * and might have to sum it in parallel runs.
1205 if (!DOMAINDECOMP(cr)
1206 || (s_min->s.ddp_count == cr->dd->ddp_count && s_b->s.ddp_count == cr->dd->ddp_count))
1208 auto fm = s_min->f.view().force();
1209 auto fb = s_b->f.view().force();
1212 /* This part of code can be incorrect with DD,
1213 * since the atom ordering in s_b and s_min might differ.
1215 for (int i = 0; i < mdatoms->homenr; i++)
1217 if (mdatoms->cFREEZE)
1219 gf = mdatoms->cFREEZE[i];
1221 for (int m = 0; m < DIM; m++)
1223 if (!opts->nFreeze[gf][m])
1225 sum += (fb[i][m] - fm[i][m]) * fb[i][m];
1232 /* We need to reorder cgs while summing */
1233 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1237 gmx_sumd(1, &sum, cr);
1240 return sum / gmx::square(s_min->fnorm);
1246 void LegacySimulator::do_cg()
1248 const char* CG = "Polak-Ribiere Conjugate Gradients";
1250 gmx_global_stat_t gstat;
1251 double tmp, minstep;
1253 real a, b, c, beta = 0.0;
1256 gmx_bool converged, foundlower;
1257 rvec mu_tot = { 0 };
1258 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1260 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1261 int m, step, nminstep;
1262 auto* mdatoms = mdAtoms->mdatoms();
1267 "Note that activating conjugate gradient energy minimization via the "
1268 "integrator .mdp option and the command gmx mdrun may "
1269 "be available in a different form in a future version of GROMACS, "
1270 "e.g. gmx minimize and an .mdp option.");
1276 // In CG, the state is extended with a search direction
1277 state_global->flags |= enumValueToBitMask(StateEntry::Cgp);
1279 // Ensure the extra per-atom state array gets allocated
1280 state_change_natoms(state_global, state_global->natoms);
1282 // Initialize the search direction to zero
1283 for (RVec& cg_p : state_global->cg_p)
1289 /* Create 4 states on the stack and extract pointers that we will swap */
1290 em_state_t s0{}, s1{}, s2{}, s3{};
1291 em_state_t* s_min = &s0;
1292 em_state_t* s_a = &s1;
1293 em_state_t* s_b = &s2;
1294 em_state_t* s_c = &s3;
1296 ObservablesReducer observablesReducer = observablesReducerBuilder->build();
1298 /* Init em and store the local state in s_min */
1317 const bool simulationsShareState = false;
1318 gmx_mdoutf* outf = init_mdoutf(fplog,
1329 StartingBehavior::NewSimulation,
1330 simulationsShareState,
1332 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
1338 StartingBehavior::NewSimulation,
1339 simulationsShareState,
1340 mdModulesNotifiers);
1342 /* Print to log file */
1343 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1345 /* Max number of steps */
1346 number_steps = inputrec->nsteps;
1350 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1354 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1357 EnergyEvaluator energyEvaluator{ fplog,
1369 &observablesReducer,
1378 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1379 /* do_force always puts the charge groups in the box and shifts again
1380 * We do not unshift, so molecules are always whole in congrad.c
1382 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE, step);
1383 observablesReducer.markAsReadyToReduce();
1387 /* Copy stuff to the energy bin for easy printing etc. */
1388 matrix nullBox = {};
1389 energyOutput.addDataAtEnergyStep(false,
1391 static_cast<double>(step),
1405 EnergyOutput::printHeader(fplog, step, step);
1406 energyOutput.printStepToEnergyFile(
1407 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
1410 /* Estimate/guess the initial stepsize */
1411 stepsize = inputrec->em_stepsize / s_min->fnorm;
1415 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1416 fprintf(stderr, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1417 fprintf(stderr, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1418 fprintf(stderr, "\n");
1419 /* and copy to the log file too... */
1420 fprintf(fplog, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1421 fprintf(fplog, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1422 fprintf(fplog, "\n");
1424 /* Start the loop over CG steps.
1425 * Each successful step is counted, and we continue until
1426 * we either converge or reach the max number of steps.
1429 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1432 /* start taking steps in a new direction
1433 * First time we enter the routine, beta=0, and the direction is
1434 * simply the negative gradient.
1437 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1438 gmx::ArrayRef<gmx::RVec> pm = s_min->s.cg_p;
1439 gmx::ArrayRef<const gmx::RVec> sfm = s_min->f.view().force();
1442 for (int i = 0; i < mdatoms->homenr; i++)
1444 if (mdatoms->cFREEZE)
1446 gf = mdatoms->cFREEZE[i];
1448 for (m = 0; m < DIM; m++)
1450 if (!inputrec->opts.nFreeze[gf][m])
1452 pm[i][m] = sfm[i][m] + beta * pm[i][m];
1453 gpa -= pm[i][m] * sfm[i][m];
1454 /* f is negative gradient, thus the sign */
1463 /* Sum the gradient along the line across CPUs */
1466 gmx_sumd(1, &gpa, cr);
1469 /* Calculate the norm of the search vector */
1470 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1472 /* Just in case stepsize reaches zero due to numerical precision... */
1475 stepsize = inputrec->em_stepsize / pnorm;
1479 * Double check the value of the derivative in the search direction.
1480 * If it is positive it must be due to the old information in the
1481 * CG formula, so just remove that and start over with beta=0.
1482 * This corresponds to a steepest descent step.
1487 step--; /* Don't count this step since we are restarting */
1488 continue; /* Go back to the beginning of the big for-loop */
1491 /* Calculate minimum allowed stepsize, before the average (norm)
1492 * relative change in coordinate is smaller than precision
1495 auto s_min_x = makeArrayRef(s_min->s.x);
1496 for (int i = 0; i < mdatoms->homenr; i++)
1498 for (m = 0; m < DIM; m++)
1500 tmp = fabs(s_min_x[i][m]);
1505 tmp = pm[i][m] / tmp;
1506 minstep += tmp * tmp;
1509 /* Add up from all CPUs */
1512 gmx_sumd(1, &minstep, cr);
1515 minstep = GMX_REAL_EPS / sqrt(minstep / (3 * top_global.natoms));
1517 if (stepsize < minstep)
1523 /* Write coordinates if necessary */
1524 do_x = do_per_step(step, inputrec->nstxout);
1525 do_f = do_per_step(step, inputrec->nstfout);
1528 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, step, s_min, state_global, observablesHistory);
1530 /* Take a step downhill.
1531 * In theory, we should minimize the function along this direction.
1532 * That is quite possible, but it turns out to take 5-10 function evaluations
1533 * for each line. However, we dont really need to find the exact minimum -
1534 * it is much better to start a new CG step in a modified direction as soon
1535 * as we are close to it. This will save a lot of energy evaluations.
1537 * In practice, we just try to take a single step.
1538 * If it worked (i.e. lowered the energy), we increase the stepsize but
1539 * the continue straight to the next CG step without trying to find any minimum.
1540 * If it didn't work (higher energy), there must be a minimum somewhere between
1541 * the old position and the new one.
1543 * Due to the finite numerical accuracy, it turns out that it is a good idea
1544 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1545 * This leads to lower final energies in the tests I've done. / Erik
1547 s_a->epot = s_min->epot;
1549 c = a + stepsize; /* reference position along line is zero */
1551 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1553 em_dd_partition_system(fplog,
1571 /* Take a trial step (new coords in s_c) */
1572 do_em_step(cr, inputrec, mdatoms, s_min, c, s_min->s.cg_p.constArrayRefWithPadding(), s_c, constr, -1);
1575 /* Calculate energy for the trial step */
1576 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE, step);
1577 observablesReducer.markAsReadyToReduce();
1579 /* Calc derivative along line */
1580 const rvec* pc = s_c->s.cg_p.rvec_array();
1581 gmx::ArrayRef<const gmx::RVec> sfc = s_c->f.view().force();
1583 for (int i = 0; i < mdatoms->homenr; i++)
1585 for (m = 0; m < DIM; m++)
1587 gpc -= pc[i][m] * sfc[i][m]; /* f is negative gradient, thus the sign */
1590 /* Sum the gradient along the line across CPUs */
1593 gmx_sumd(1, &gpc, cr);
1596 /* This is the max amount of increase in energy we tolerate */
1597 tmp = std::sqrt(GMX_REAL_EPS) * fabs(s_a->epot);
1599 /* Accept the step if the energy is lower, or if it is not significantly higher
1600 * and the line derivative is still negative.
1602 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1605 /* Great, we found a better energy. Increase step for next iteration
1606 * if we are still going down, decrease it otherwise
1610 stepsize *= 1.618034; /* The golden section */
1614 stepsize *= 0.618034; /* 1/golden section */
1619 /* New energy is the same or higher. We will have to do some work
1620 * to find a smaller value in the interval. Take smaller step next time!
1623 stepsize *= 0.618034;
1627 /* OK, if we didn't find a lower value we will have to locate one now - there must
1628 * be one in the interval [a=0,c].
1629 * The same thing is valid here, though: Don't spend dozens of iterations to find
1630 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1631 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1633 * I also have a safeguard for potentially really pathological functions so we never
1634 * take more than 20 steps before we give up ...
1636 * If we already found a lower value we just skip this step and continue to the update.
1645 /* Select a new trial point.
1646 * If the derivatives at points a & c have different sign we interpolate to zero,
1647 * otherwise just do a bisection.
1649 if (gpa < 0 && gpc > 0)
1651 b = a + gpa * (a - c) / (gpc - gpa);
1658 /* safeguard if interpolation close to machine accuracy causes errors:
1659 * never go outside the interval
1661 if (b <= a || b >= c)
1666 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1668 /* Reload the old state */
1669 em_dd_partition_system(fplog,
1687 /* Take a trial step to this new point - new coords in s_b */
1688 do_em_step(cr, inputrec, mdatoms, s_min, b, s_min->s.cg_p.constArrayRefWithPadding(), s_b, constr, -1);
1691 /* Calculate energy for the trial step */
1692 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE, step);
1693 observablesReducer.markAsReadyToReduce();
1695 /* p does not change within a step, but since the domain decomposition
1696 * might change, we have to use cg_p of s_b here.
1698 const rvec* pb = s_b->s.cg_p.rvec_array();
1699 gmx::ArrayRef<const gmx::RVec> sfb = s_b->f.view().force();
1701 for (int i = 0; i < mdatoms->homenr; i++)
1703 for (m = 0; m < DIM; m++)
1705 gpb -= pb[i][m] * sfb[i][m]; /* f is negative gradient, thus the sign */
1708 /* Sum the gradient along the line across CPUs */
1711 gmx_sumd(1, &gpb, cr);
1716 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a->epot, s_b->epot, s_c->epot, gpb);
1719 epot_repl = s_b->epot;
1721 /* Keep one of the intervals based on the value of the derivative at the new point */
1724 /* Replace c endpoint with b */
1725 swap_em_state(&s_b, &s_c);
1731 /* Replace a endpoint with b */
1732 swap_em_state(&s_b, &s_a);
1738 * Stop search as soon as we find a value smaller than the endpoints.
1739 * Never run more than 20 steps, no matter what.
1742 } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) && (nminstep < 20));
1744 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot) * GMX_REAL_EPS || nminstep >= 20)
1746 /* OK. We couldn't find a significantly lower energy.
1747 * If beta==0 this was steepest descent, and then we give up.
1748 * If not, set beta=0 and restart with steepest descent before quitting.
1758 /* Reset memory before giving up */
1764 /* Select min energy state of A & C, put the best in B.
1766 if (s_c->epot < s_a->epot)
1770 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c->epot, s_a->epot);
1772 swap_em_state(&s_b, &s_c);
1779 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a->epot, s_c->epot);
1781 swap_em_state(&s_b, &s_a);
1789 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n", s_c->epot);
1791 swap_em_state(&s_b, &s_c);
1795 /* new search direction */
1796 /* beta = 0 means forget all memory and restart with steepest descents. */
1797 if (nstcg && ((step % nstcg) == 0))
1803 /* s_min->fnorm cannot be zero, because then we would have converged
1807 /* Polak-Ribiere update.
1808 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1810 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1812 /* Limit beta to prevent oscillations */
1813 if (fabs(beta) > 5.0)
1819 /* update positions */
1820 swap_em_state(&s_min, &s_b);
1823 /* Print it if necessary */
1826 if (mdrunOptions.verbose)
1828 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1830 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1833 s_min->fnorm / sqrtNumAtoms,
1838 /* Store the new (lower) energies */
1839 matrix nullBox = {};
1840 energyOutput.addDataAtEnergyStep(false,
1842 static_cast<double>(step),
1856 do_log = do_per_step(step, inputrec->nstlog);
1857 do_ene = do_per_step(step, inputrec->nstenergy);
1859 imdSession->fillEnergyRecord(step, TRUE);
1863 EnergyOutput::printHeader(fplog, step, step);
1865 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1869 do_log ? fplog : nullptr,
1876 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1877 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x, 0))
1879 imdSession->sendPositionsAndEnergies();
1882 /* Stop when the maximum force lies below tolerance.
1883 * If we have reached machine precision, converged is already set to true.
1885 converged = converged || (s_min->fmax < inputrec->em_tol);
1886 observablesReducer.markAsReadyToReduce();
1887 } /* End of the loop */
1891 step--; /* we never took that last step in this case */
1893 if (s_min->fmax > inputrec->em_tol)
1897 warn_step(fplog, inputrec->em_tol, s_min->fmax, step - 1 == number_steps, FALSE);
1904 /* If we printed energy and/or logfile last step (which was the last step)
1905 * we don't have to do it again, but otherwise print the final values.
1909 /* Write final value to log since we didn't do anything the last step */
1910 EnergyOutput::printHeader(fplog, step, step);
1912 if (!do_ene || !do_log)
1914 /* Write final energy file entries */
1915 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1919 !do_log ? fplog : nullptr,
1927 /* Print some stuff... */
1930 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1934 * For accurate normal mode calculation it is imperative that we
1935 * store the last conformation into the full precision binary trajectory.
1937 * However, we should only do it if we did NOT already write this step
1938 * above (which we did if do_x or do_f was true).
1940 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1941 * in the trajectory with the same step number.
1943 do_x = !do_per_step(step, inputrec->nstxout);
1944 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1947 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, s_min, state_global, observablesHistory);
1952 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1953 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1954 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1956 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1959 finish_em(cr, outf, walltime_accounting, wcycle);
1961 /* To print the actual number of steps we needed somewhere */
1962 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1966 void LegacySimulator::do_lbfgs()
1968 static const char* LBFGS = "Low-Memory BFGS Minimizer";
1970 gmx_global_stat_t gstat;
1971 auto* mdatoms = mdAtoms->mdatoms();
1976 "Note that activating L-BFGS energy minimization via the "
1977 "integrator .mdp option and the command gmx mdrun may "
1978 "be available in a different form in a future version of GROMACS, "
1979 "e.g. gmx minimize and an .mdp option.");
1981 if (DOMAINDECOMP(cr))
1983 gmx_fatal(FARGS, "L_BFGS is currently not supported");
1987 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1990 if (nullptr != constr)
1994 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1995 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1998 const int n = 3 * state_global->natoms;
1999 const int nmaxcorr = inputrec->nbfgscorr;
2001 std::vector<real> p(n);
2002 std::vector<real> rho(nmaxcorr);
2003 std::vector<real> alpha(nmaxcorr);
2005 std::vector<std::vector<real>> dx(nmaxcorr);
2006 for (auto& dxCorr : dx)
2011 std::vector<std::vector<real>> dg(nmaxcorr);
2012 for (auto& dgCorr : dg)
2020 ObservablesReducer observablesReducer = observablesReducerBuilder->build();
2041 const bool simulationsShareState = false;
2042 gmx_mdoutf* outf = init_mdoutf(fplog,
2053 StartingBehavior::NewSimulation,
2054 simulationsShareState,
2056 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
2062 StartingBehavior::NewSimulation,
2063 simulationsShareState,
2064 mdModulesNotifiers);
2066 const int start = 0;
2067 const int end = mdatoms->homenr;
2069 /* We need 4 working states */
2070 em_state_t s0{}, s1{}, s2{}, s3{};
2071 em_state_t* sa = &s0;
2072 em_state_t* sb = &s1;
2073 em_state_t* sc = &s2;
2074 em_state_t* last = &s3;
2075 /* Initialize by copying the state from ems (we could skip x and f here) */
2080 /* Print to log file */
2081 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
2083 /* Max number of steps */
2084 const int number_steps = inputrec->nsteps;
2086 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
2087 std::vector<bool> frozen(n);
2089 for (int i = start; i < end; i++)
2091 if (mdatoms->cFREEZE)
2093 gf = mdatoms->cFREEZE[i];
2095 for (int m = 0; m < DIM; m++)
2097 frozen[3 * i + m] = (inputrec->opts.nFreeze[gf][m] != 0);
2102 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
2106 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
2111 vsite->construct(state_global->x, {}, state_global->box, VSiteOperation::Positions);
2114 /* Call the force routine and some auxiliary (neighboursearching etc.) */
2115 /* do_force always puts the charge groups in the box and shifts again
2116 * We do not unshift, so molecules are always whole
2119 EnergyEvaluator energyEvaluator{ fplog,
2131 &observablesReducer,
2141 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE, step);
2145 /* Copy stuff to the energy bin for easy printing etc. */
2146 matrix nullBox = {};
2147 energyOutput.addDataAtEnergyStep(false,
2149 static_cast<double>(step),
2163 EnergyOutput::printHeader(fplog, step, step);
2164 energyOutput.printStepToEnergyFile(
2165 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
2168 /* Set the initial step.
2169 * since it will be multiplied by the non-normalized search direction
2170 * vector (force vector the first time), we scale it by the
2171 * norm of the force.
2176 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2177 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2178 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2179 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2180 fprintf(stderr, "\n");
2181 /* and copy to the log file too... */
2182 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2183 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2184 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2185 fprintf(fplog, "\n");
2188 // Point is an index to the memory of search directions, where 0 is the first one.
2191 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
2192 real* fInit = static_cast<real*>(ems.f.view().force().data()[0]);
2193 for (int i = 0; i < n; i++)
2197 dx[point][i] = fInit[i]; /* Initial search direction */
2205 // Stepsize will be modified during the search, and actually it is not critical
2206 // (the main efficiency in the algorithm comes from changing directions), but
2207 // we still need an initial value, so estimate it as the inverse of the norm
2208 // so we take small steps where the potential fluctuates a lot.
2209 double stepsize = 1.0 / ems.fnorm;
2211 /* Start the loop over BFGS steps.
2212 * Each successful step is counted, and we continue until
2213 * we either converge or reach the max number of steps.
2221 /* Set the gradient from the force */
2222 bool converged = false;
2223 for (int step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
2226 /* Write coordinates if necessary */
2227 const bool do_x = do_per_step(step, inputrec->nstxout);
2228 const bool do_f = do_per_step(step, inputrec->nstfout);
2233 mdof_flags |= MDOF_X;
2238 mdof_flags |= MDOF_F;
2243 mdof_flags |= MDOF_IMD;
2246 gmx::WriteCheckpointDataHolder checkpointDataHolder;
2247 mdoutf_write_to_trajectory_files(fplog,
2253 static_cast<real>(step),
2257 ems.f.view().force(),
2258 &checkpointDataHolder);
2260 /* Do the linesearching in the direction dx[point][0..(n-1)] */
2262 /* make s a pointer to current search direction - point=0 first time we get here */
2263 gmx::ArrayRef<const real> s = dx[point];
2265 const real* xx = static_cast<real*>(ems.s.x.rvec_array()[0]);
2266 const real* ff = static_cast<real*>(ems.f.view().force().data()[0]);
2268 // calculate line gradient in position A
2270 for (int i = 0; i < n; i++)
2272 gpa -= s[i] * ff[i];
2275 /* Calculate minimum allowed stepsize along the line, before the average (norm)
2276 * relative change in coordinate is smaller than precision
2279 for (int i = 0; i < n; i++)
2281 double tmp = fabs(xx[i]);
2287 minstep += tmp * tmp;
2289 minstep = GMX_REAL_EPS / sqrt(minstep / n);
2291 if (stepsize < minstep)
2297 // Before taking any steps along the line, store the old position
2299 real* lastx = static_cast<real*>(last->s.x.data()[0]);
2300 real* lastf = static_cast<real*>(last->f.view().force().data()[0]);
2301 const real Epot0 = ems.epot;
2305 /* Take a step downhill.
2306 * In theory, we should find the actual minimum of the function in this
2307 * direction, somewhere along the line.
2308 * That is quite possible, but it turns out to take 5-10 function evaluations
2309 * for each line. However, we dont really need to find the exact minimum -
2310 * it is much better to start a new BFGS step in a modified direction as soon
2311 * as we are close to it. This will save a lot of energy evaluations.
2313 * In practice, we just try to take a single step.
2314 * If it worked (i.e. lowered the energy), we increase the stepsize but
2315 * continue straight to the next BFGS step without trying to find any minimum,
2316 * i.e. we change the search direction too. If the line was smooth, it is
2317 * likely we are in a smooth region, and then it makes sense to take longer
2318 * steps in the modified search direction too.
2320 * If it didn't work (higher energy), there must be a minimum somewhere between
2321 * the old position and the new one. Then we need to start by finding a lower
2322 * value before we change search direction. Since the energy was apparently
2323 * quite rough, we need to decrease the step size.
2325 * Due to the finite numerical accuracy, it turns out that it is a good idea
2326 * to accept a SMALL increase in energy, if the derivative is still downhill.
2327 * This leads to lower final energies in the tests I've done. / Erik
2330 // State "A" is the first position along the line.
2331 // reference position along line is initially zero
2334 // Check stepsize first. We do not allow displacements
2335 // larger than emstep.
2341 // Pick a new position C by adding stepsize to A.
2344 // Calculate what the largest change in any individual coordinate
2345 // would be (translation along line * gradient along line)
2347 for (int i = 0; i < n; i++)
2349 real delta = c * s[i];
2350 if (delta > maxdelta)
2355 // If any displacement is larger than the stepsize limit, reduce the step
2356 if (maxdelta > inputrec->em_stepsize)
2360 } while (maxdelta > inputrec->em_stepsize);
2362 // Take a trial step and move the coordinate array xc[] to position C
2363 real* xc = static_cast<real*>(sc->s.x.rvec_array()[0]);
2364 for (int i = 0; i < n; i++)
2366 xc[i] = lastx[i] + c * s[i];
2370 // Calculate energy for the trial step in position C
2371 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE, step);
2373 // Calc line gradient in position C
2374 real* fc = static_cast<real*>(sc->f.view().force()[0]);
2376 for (int i = 0; i < n; i++)
2378 gpc -= s[i] * fc[i]; /* f is negative gradient, thus the sign */
2380 /* Sum the gradient along the line across CPUs */
2383 gmx_sumd(1, &gpc, cr);
2386 // This is the max amount of increase in energy we tolerate.
2387 // By allowing VERY small changes (close to numerical precision) we
2388 // frequently find even better (lower) final energies.
2389 double tmp = std::sqrt(GMX_REAL_EPS) * fabs(sa->epot);
2391 // Accept the step if the energy is lower in the new position C (compared to A),
2392 // or if it is not significantly higher and the line derivative is still negative.
2393 bool foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2394 // If true, great, we found a better energy. We no longer try to alter the
2395 // stepsize, but simply accept this new better position. The we select a new
2396 // search direction instead, which will be much more efficient than continuing
2397 // to take smaller steps along a line. Set fnorm based on the new C position,
2398 // which will be used to update the stepsize to 1/fnorm further down.
2400 // If false, the energy is NOT lower in point C, i.e. it will be the same
2401 // or higher than in point A. In this case it is pointless to move to point C,
2402 // so we will have to do more iterations along the same line to find a smaller
2403 // value in the interval [A=0.0,C].
2404 // Here, A is still 0.0, but that will change when we do a search in the interval
2405 // [0.0,C] below. That search we will do by interpolation or bisection rather
2406 // than with the stepsize, so no need to modify it. For the next search direction
2407 // it will be reset to 1/fnorm anyway.
2412 // OK, if we didn't find a lower value we will have to locate one now - there must
2413 // be one in the interval [a,c].
2414 // The same thing is valid here, though: Don't spend dozens of iterations to find
2415 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2416 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2417 // I also have a safeguard for potentially really pathological functions so we never
2418 // take more than 20 steps before we give up.
2419 // If we already found a lower value we just skip this step and continue to the update.
2424 // Select a new trial point B in the interval [A,C].
2425 // If the derivatives at points a & c have different sign we interpolate to zero,
2426 // otherwise just do a bisection since there might be multiple minima/maxima
2427 // inside the interval.
2429 if (gpa < 0 && gpc > 0)
2431 b = a + gpa * (a - c) / (gpc - gpa);
2438 /* safeguard if interpolation close to machine accuracy causes errors:
2439 * never go outside the interval
2441 if (b <= a || b >= c)
2446 // Take a trial step to point B
2447 real* xb = static_cast<real*>(sb->s.x.rvec_array()[0]);
2448 for (int i = 0; i < n; i++)
2450 xb[i] = lastx[i] + b * s[i];
2454 // Calculate energy for the trial step in point B
2455 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE, step);
2458 // Calculate gradient in point B
2459 real* fb = static_cast<real*>(sb->f.view().force()[0]);
2461 for (int i = 0; i < n; i++)
2463 gpb -= s[i] * fb[i]; /* f is negative gradient, thus the sign */
2465 /* Sum the gradient along the line across CPUs */
2468 gmx_sumd(1, &gpb, cr);
2471 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2472 // at the new point B, and rename the endpoints of this new interval A and C.
2475 /* Replace c endpoint with b */
2477 /* copy state b to c */
2482 /* Replace a endpoint with b */
2484 /* copy state b to a */
2489 * Stop search as soon as we find a value smaller than the endpoints,
2490 * or if the tolerance is below machine precision.
2491 * Never run more than 20 steps, no matter what.
2494 } while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2496 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2498 /* OK. We couldn't find a significantly lower energy.
2499 * If ncorr==0 this was steepest descent, and then we give up.
2500 * If not, reset memory to restart as steepest descent before quitting.
2512 /* Search in gradient direction */
2513 for (int i = 0; i < n; i++)
2515 dx[point][i] = ff[i];
2517 /* Reset stepsize */
2518 stepsize = 1.0 / fnorm;
2523 /* Select min energy state of A & C, put the best in xx/ff/Epot
2525 if (sc->epot < sa->epot)
2546 /* Update the memory information, and calculate a new
2547 * approximation of the inverse hessian
2550 /* Have new data in Epot, xx, ff */
2551 if (ncorr < nmaxcorr)
2556 for (int i = 0; i < n; i++)
2558 dg[point][i] = lastf[i] - ff[i];
2559 dx[point][i] *= step_taken;
2564 for (int i = 0; i < n; i++)
2566 dgdg += dg[point][i] * dg[point][i];
2567 dgdx += dg[point][i] * dx[point][i];
2570 const real diag = dgdx / dgdg;
2572 rho[point] = 1.0 / dgdx;
2575 if (point >= nmaxcorr)
2581 for (int i = 0; i < n; i++)
2588 /* Recursive update. First go back over the memory points */
2589 for (int k = 0; k < ncorr; k++)
2598 for (int i = 0; i < n; i++)
2600 sq += dx[cp][i] * p[i];
2603 alpha[cp] = rho[cp] * sq;
2605 for (int i = 0; i < n; i++)
2607 p[i] -= alpha[cp] * dg[cp][i];
2611 for (int i = 0; i < n; i++)
2616 /* And then go forward again */
2617 for (int k = 0; k < ncorr; k++)
2620 for (int i = 0; i < n; i++)
2622 yr += p[i] * dg[cp][i];
2625 real beta = rho[cp] * yr;
2626 beta = alpha[cp] - beta;
2628 for (int i = 0; i < n; i++)
2630 p[i] += beta * dx[cp][i];
2640 for (int i = 0; i < n; i++)
2644 dx[point][i] = p[i];
2652 /* Print it if necessary */
2655 if (mdrunOptions.verbose)
2657 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2659 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2662 ems.fnorm / sqrtNumAtoms,
2667 /* Store the new (lower) energies */
2668 matrix nullBox = {};
2669 energyOutput.addDataAtEnergyStep(false,
2671 static_cast<double>(step),
2685 do_log = do_per_step(step, inputrec->nstlog);
2686 do_ene = do_per_step(step, inputrec->nstenergy);
2688 imdSession->fillEnergyRecord(step, TRUE);
2692 EnergyOutput::printHeader(fplog, step, step);
2694 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2698 do_log ? fplog : nullptr,
2705 /* Send x and E to IMD client, if bIMD is TRUE. */
2706 if (imdSession->run(step, TRUE, state_global->box, state_global->x, 0) && MASTER(cr))
2708 imdSession->sendPositionsAndEnergies();
2711 // Reset stepsize in we are doing more iterations
2714 /* Stop when the maximum force lies below tolerance.
2715 * If we have reached machine precision, converged is already set to true.
2717 converged = converged || (ems.fmax < inputrec->em_tol);
2718 observablesReducer.markAsReadyToReduce();
2719 } /* End of the loop */
2723 step--; /* we never took that last step in this case */
2725 if (ems.fmax > inputrec->em_tol)
2729 warn_step(fplog, inputrec->em_tol, ems.fmax, step - 1 == number_steps, FALSE);
2734 /* If we printed energy and/or logfile last step (which was the last step)
2735 * we don't have to do it again, but otherwise print the final values.
2737 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2739 EnergyOutput::printHeader(fplog, step, step);
2741 if (!do_ene || !do_log) /* Write final energy file entries */
2743 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2747 !do_log ? fplog : nullptr,
2754 /* Print some stuff... */
2757 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2761 * For accurate normal mode calculation it is imperative that we
2762 * store the last conformation into the full precision binary trajectory.
2764 * However, we should only do it if we did NOT already write this step
2765 * above (which we did if do_x or do_f was true).
2767 const bool do_x = !do_per_step(step, inputrec->nstxout);
2768 const bool do_f = !do_per_step(step, inputrec->nstfout);
2770 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, &ems, state_global, observablesHistory);
2774 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2775 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2776 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2778 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2781 finish_em(cr, outf, walltime_accounting, wcycle);
2783 /* To print the actual number of steps we needed somewhere */
2784 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2787 void LegacySimulator::do_steep()
2789 const char* SD = "Steepest Descents";
2790 gmx_global_stat_t gstat;
2793 gmx_bool bDone, bAbort, do_x, do_f;
2795 rvec mu_tot = { 0 };
2798 int steps_accepted = 0;
2799 auto* mdatoms = mdAtoms->mdatoms();
2804 "Note that activating steepest-descent energy minimization via the "
2805 "integrator .mdp option and the command gmx mdrun may "
2806 "be available in a different form in a future version of GROMACS, "
2807 "e.g. gmx minimize and an .mdp option.");
2809 /* Create 2 states on the stack and extract pointers that we will swap */
2810 em_state_t s0{}, s1{};
2811 em_state_t* s_min = &s0;
2812 em_state_t* s_try = &s1;
2814 ObservablesReducer observablesReducer = observablesReducerBuilder->build();
2816 /* Init em and store the local state in s_try */
2835 const bool simulationsShareState = false;
2836 gmx_mdoutf* outf = init_mdoutf(fplog,
2847 StartingBehavior::NewSimulation,
2848 simulationsShareState,
2850 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
2856 StartingBehavior::NewSimulation,
2857 simulationsShareState,
2858 mdModulesNotifiers);
2860 /* Print to log file */
2861 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2863 /* Set variables for stepsize (in nm). This is the largest
2864 * step that we are going to make in any direction.
2866 ustep = inputrec->em_stepsize;
2869 /* Max number of steps */
2870 nsteps = inputrec->nsteps;
2874 /* Print to the screen */
2875 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2879 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2881 EnergyEvaluator energyEvaluator{ fplog,
2893 &observablesReducer,
2901 /**** HERE STARTS THE LOOP ****
2902 * count is the counter for the number of steps
2903 * bDone will be TRUE when the minimization has converged
2904 * bAbort will be TRUE when nsteps steps have been performed or when
2905 * the stepsize becomes smaller than is reasonable for machine precision
2910 while (!bDone && !bAbort)
2912 bAbort = (nsteps >= 0) && (count == nsteps);
2914 /* set new coordinates, except for first step */
2915 bool validStep = true;
2918 validStep = do_em_step(
2919 cr, inputrec, mdatoms, s_min, stepsize, s_min->f.view().forceWithPadding(), s_try, constr, count);
2924 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0, count);
2928 // Signal constraint error during stepping with energy=inf
2929 s_try->epot = std::numeric_limits<real>::infinity();
2934 EnergyOutput::printHeader(fplog, count, count);
2939 s_min->epot = s_try->epot;
2942 /* Print it if necessary */
2945 if (mdrunOptions.verbose)
2948 "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2954 ((count == 0) || (s_try->epot < s_min->epot)) ? '\n' : '\r');
2958 if ((count == 0) || (s_try->epot < s_min->epot))
2960 /* Store the new (lower) energies */
2961 matrix nullBox = {};
2962 energyOutput.addDataAtEnergyStep(false,
2964 static_cast<double>(count),
2978 imdSession->fillEnergyRecord(count, TRUE);
2980 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2981 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2982 energyOutput.printStepToEnergyFile(
2983 mdoutf_get_fp_ene(outf), TRUE, do_dr, do_or, fplog, count, count, fr->fcdata.get(), nullptr);
2988 /* Now if the new energy is smaller than the previous...
2989 * or if this is the first step!
2990 * or if we did random steps!
2993 if ((count == 0) || (s_try->epot < s_min->epot))
2997 /* Test whether the convergence criterion is met... */
2998 bDone = (s_try->fmax < inputrec->em_tol);
3000 /* Copy the arrays for force, positions and energy */
3001 /* The 'Min' array always holds the coords and forces of the minimal
3003 swap_em_state(&s_min, &s_try);
3009 /* Write to trn, if necessary */
3010 do_x = do_per_step(steps_accepted, inputrec->nstxout);
3011 do_f = do_per_step(steps_accepted, inputrec->nstfout);
3013 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, count, s_min, state_global, observablesHistory);
3017 /* If energy is not smaller make the step smaller... */
3020 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
3022 /* Reload the old state */
3023 em_dd_partition_system(fplog,
3042 // If the force is very small after finishing minimization,
3043 // we risk dividing by zero when calculating the step size.
3044 // So we check first if the minimization has stopped before
3045 // trying to obtain a new step size.
3048 /* Determine new step */
3049 stepsize = ustep / s_min->fmax;
3052 /* Check if stepsize is too small, with 1 nm as a characteristic length */
3054 if (count == nsteps || ustep < 1e-12)
3056 if (count == nsteps || ustep < 1e-6)
3061 warn_step(fplog, inputrec->em_tol, s_min->fmax, count == nsteps, constr != nullptr);
3066 /* Send IMD energies and positions, if bIMD is TRUE. */
3067 if (imdSession->run(count,
3069 MASTER(cr) ? state_global->box : nullptr,
3070 MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>(),
3074 imdSession->sendPositionsAndEnergies();
3078 observablesReducer.markAsReadyToReduce();
3079 } /* End of the loop */
3081 /* Print some data... */
3084 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
3086 write_em_traj(fplog,
3090 inputrec->nstfout != 0,
3091 ftp2fn(efSTO, nfile, fnm),
3097 observablesHistory);
3101 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
3103 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
3104 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
3107 finish_em(cr, outf, walltime_accounting, wcycle);
3109 /* To print the actual number of steps we needed somewhere */
3111 // TODO: Avoid changing inputrec (#3854)
3112 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
3113 nonConstInputrec->nsteps = count;
3116 walltime_accounting_set_nsteps_done(walltime_accounting, count);
3119 void LegacySimulator::do_nm()
3121 const char* NM = "Normal Mode Analysis";
3123 gmx_global_stat_t gstat;
3125 rvec mu_tot = { 0 };
3127 gmx_bool bSparse; /* use sparse matrix storage format */
3129 gmx_sparsematrix_t* sparse_matrix = nullptr;
3130 real* full_matrix = nullptr;
3132 /* added with respect to mdrun */
3134 real der_range = 10.0 * std::sqrt(GMX_REAL_EPS);
3136 bool bIsMaster = MASTER(cr);
3137 auto* mdatoms = mdAtoms->mdatoms();
3142 "Note that activating normal-mode analysis via the integrator "
3143 ".mdp option and the command gmx mdrun may "
3144 "be available in a different form in a future version of GROMACS, "
3145 "e.g. gmx normal-modes.");
3147 if (constr != nullptr)
3151 "Constraints present with Normal Mode Analysis, this combination is not supported");
3154 gmx_shellfc_t* shellfc;
3156 em_state_t state_work{};
3158 fr->longRangeNonbondeds->updateAfterPartition(*mdAtoms->mdatoms());
3159 ObservablesReducer observablesReducer = observablesReducerBuilder->build();
3161 /* Init em and store the local state in state_minimum */
3180 const bool simulationsShareState = false;
3181 gmx_mdoutf* outf = init_mdoutf(fplog,
3192 StartingBehavior::NewSimulation,
3193 simulationsShareState,
3196 std::vector<int> atom_index = get_atom_index(top_global);
3197 std::vector<gmx::RVec> fneg(atom_index.size(), { 0, 0, 0 });
3198 snew(dfdx, atom_index.size());
3204 "NOTE: This version of GROMACS has been compiled in single precision,\n"
3205 " which MIGHT not be accurate enough for normal mode analysis.\n"
3206 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
3207 " are fairly modest even if you recompile in double precision.\n\n");
3211 /* Check if we can/should use sparse storage format.
3213 * Sparse format is only useful when the Hessian itself is sparse, which it
3214 * will be when we use a cutoff.
3215 * For small systems (n<1000) it is easier to always use full matrix format, though.
3217 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
3219 GMX_LOG(mdlog.warning)
3220 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
3223 else if (atom_index.size() < 1000)
3225 GMX_LOG(mdlog.warning)
3226 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
3232 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
3236 /* Number of dimensions, based on real atoms, that is not vsites or shell */
3237 sz = DIM * atom_index.size();
3239 fprintf(stderr, "Allocating Hessian memory...\n\n");
3243 sparse_matrix = gmx_sparsematrix_init(sz);
3244 sparse_matrix->compressed_symmetric = TRUE;
3248 snew(full_matrix, sz * sz);
3251 /* Write start time and temperature */
3252 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
3254 /* fudge nr of steps to nr of atoms */
3256 // TODO: Avoid changing inputrec (#3854)
3257 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
3258 nonConstInputrec->nsteps = atom_index.size() * 2;
3264 "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
3269 nnodes = cr->nnodes;
3271 /* Make evaluate_energy do a single node force calculation */
3273 EnergyEvaluator energyEvaluator{ fplog,
3285 &observablesReducer,
3292 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE, 0);
3293 cr->nnodes = nnodes;
3295 /* if forces are not small, warn user */
3296 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
3298 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
3299 if (state_work.fmax > 1.0e-3)
3301 GMX_LOG(mdlog.warning)
3303 "The force is probably not small enough to "
3304 "ensure that you are at a minimum.\n"
3305 "Be aware that negative eigenvalues may occur\n"
3306 "when the resulting matrix is diagonalized.");
3309 /***********************************************************
3311 * Loop over all pairs in matrix
3313 * do_force called twice. Once with positive and
3314 * once with negative displacement
3316 ************************************************************/
3318 /* Steps are divided one by one over the nodes */
3320 auto state_work_x = makeArrayRef(state_work.s.x);
3321 auto state_work_f = state_work.f.view().force();
3322 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
3324 size_t atom = atom_index[aid];
3325 for (size_t d = 0; d < DIM; d++)
3328 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
3331 x_min = state_work_x[atom][d];
3333 for (unsigned int dx = 0; (dx < 2); dx++)
3337 state_work_x[atom][d] = x_min - der_range;
3341 state_work_x[atom][d] = x_min + der_range;
3344 /* Make evaluate_energy do a single node force calculation */
3348 /* Now is the time to relax the shells */
3349 relax_shell_flexcon(fplog,
3352 mdrunOptions.verbose,
3363 state_work.s.natoms,
3364 state_work.s.x.arrayRefWithPadding(),
3365 state_work.s.v.arrayRefWithPadding(),
3367 state_work.s.lambda,
3369 &state_work.f.view(),
3372 fr->longRangeNonbondeds.get(),
3381 DDBalanceRegionHandler(nullptr));
3387 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid * 2 + dx, FALSE, step);
3390 cr->nnodes = nnodes;
3394 std::copy(state_work_f.begin(), state_work_f.begin() + atom_index.size(), fneg.begin());
3398 /* x is restored to original */
3399 state_work_x[atom][d] = x_min;
3401 for (size_t j = 0; j < atom_index.size(); j++)
3403 for (size_t k = 0; (k < DIM); k++)
3405 dfdx[j][k] = -(state_work_f[atom_index[j]][k] - fneg[j][k]) / (2 * der_range);
3412 # define mpi_type GMX_MPI_REAL
3413 MPI_Send(dfdx[0], atom_index.size() * DIM, mpi_type, MASTER(cr), cr->nodeid, cr->mpi_comm_mygroup);
3418 for (index node = 0; (node < nnodes && aid + node < ssize(atom_index)); node++)
3424 MPI_Recv(dfdx[0], atom_index.size() * DIM, mpi_type, node, node, cr->mpi_comm_mygroup, &stat);
3429 row = (aid + node) * DIM + d;
3431 for (size_t j = 0; j < atom_index.size(); j++)
3433 for (size_t k = 0; k < DIM; k++)
3439 if (col >= row && dfdx[j][k] != 0.0)
3441 gmx_sparsematrix_increment_value(sparse_matrix, row, col, dfdx[j][k]);
3446 full_matrix[row * sz + col] = dfdx[j][k];
3453 if (mdrunOptions.verbose && fplog)
3458 /* write progress */
3459 if (bIsMaster && mdrunOptions.verbose)
3462 "\rFinished step %d out of %td",
3463 std::min<int>(atom + nnodes, atom_index.size()),
3471 fprintf(stderr, "\n\nWriting Hessian...\n");
3472 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
3475 finish_em(cr, outf, walltime_accounting, wcycle);
3477 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size() * 2);