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>();
398 fplog, *ir, gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc), MASTER(cr), fep_state, lambda);
400 if (ir->eI == IntegrationAlgorithm::NM)
402 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
404 *shellfc = init_shell_flexcon(stdout,
406 constr ? constr->numFlexibleConstraints() : 0,
409 thisRankHasDuty(cr, DUTY_PME));
413 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI),
414 "This else currently only handles energy minimizers, consider if your algorithm "
415 "needs shell/flexible-constraint support");
417 /* With energy minimization, shells and flexible constraints are
418 * automatically minimized when treated like normal DOFS.
420 if (shellfc != nullptr)
426 if (DOMAINDECOMP(cr))
428 dd_init_local_state(*cr->dd, state_global, &ems->s);
430 /* Distribute the charge groups over the nodes from the master node */
431 dd_partition_system(fplog,
452 dd_store_state(*cr->dd, &ems->s);
456 state_change_natoms(state_global, state_global->natoms);
457 /* Just copy the state */
458 ems->s = *state_global;
459 state_change_natoms(&ems->s, ems->s.natoms);
461 mdAlgorithmsSetupAtomData(
462 cr, *ir, top_global, top, fr, &ems->f, mdAtoms, constr, vsite, shellfc ? *shellfc : nullptr);
465 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[FreeEnergyPerturbationCouplingType::Mass]);
469 // TODO how should this cross-module support dependency be managed?
470 if (ir->eConstrAlg == ConstraintAlgorithm::Shake && gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
473 "Can not do energy minimization with %s, use %s\n",
474 enumValueToString(ConstraintAlgorithm::Shake),
475 enumValueToString(ConstraintAlgorithm::Lincs));
478 if (!ir->bContinuation)
480 /* Constrain the starting coordinates */
481 bool needsLogging = true;
482 bool computeEnergy = true;
483 bool computeVirial = false;
485 constr->apply(needsLogging,
490 ems->s.x.arrayRefWithPadding(),
491 ems->s.x.arrayRefWithPadding(),
494 ems->s.lambda[FreeEnergyPerturbationCouplingType::Fep],
496 gmx::ArrayRefWithPadding<RVec>(),
499 gmx::ConstraintVariable::Positions);
505 *gstat = global_stat_init(ir);
512 calc_shifts(ems->s.box, fr->shift_vec);
515 //! Finalize the minimization
516 static void finish_em(const t_commrec* cr,
518 gmx_walltime_accounting_t walltime_accounting,
519 gmx_wallcycle* wcycle)
521 if (!thisRankHasDuty(cr, DUTY_PME))
523 /* Tell the PME only node to finish */
524 gmx_pme_send_finish(cr);
529 em_time_end(walltime_accounting, wcycle);
532 //! Swap two different EM states during minimization
533 static void swap_em_state(em_state_t** ems1, em_state_t** ems2)
542 //! Save the EM trajectory
543 static void write_em_traj(FILE* fplog,
549 const gmx_mtop_t& top_global,
550 const t_inputrec* ir,
553 t_state* state_global,
554 ObservablesHistory* observablesHistory)
560 mdof_flags |= MDOF_X;
564 mdof_flags |= MDOF_F;
567 /* If we want IMD output, set appropriate MDOF flag */
570 mdof_flags |= MDOF_IMD;
573 gmx::WriteCheckpointDataHolder checkpointDataHolder;
574 mdoutf_write_to_trajectory_files(fplog,
580 static_cast<double>(step),
584 state->f.view().force(),
585 &checkpointDataHolder);
587 if (confout != nullptr)
589 if (DOMAINDECOMP(cr))
591 /* If bX=true, x was collected to state_global in the call above */
594 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
596 cr->dd, state->s.ddp_count, state->s.ddp_count_cg_gl, state->s.cg_gl, state->s.x, globalXRef);
601 /* Copy the local state pointer */
602 state_global = &state->s;
607 if (ir->pbcType != PbcType::No && !ir->bPeriodicMols && DOMAINDECOMP(cr))
609 /* Make molecules whole only for confout writing */
610 do_pbc_mtop(ir->pbcType, state->s.box, &top_global, state_global->x.rvec_array());
613 write_sto_conf_mtop(confout,
616 state_global->x.rvec_array(),
624 //! \brief Do one minimization step
626 // \returns true when the step succeeded, false when a constraint error occurred
627 static bool do_em_step(const t_commrec* cr,
628 const t_inputrec* ir,
632 gmx::ArrayRefWithPadding<const gmx::RVec> force,
634 gmx::Constraints* constr,
641 int nthreads gmx_unused;
643 bool validStep = true;
648 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
650 gmx_incons("state mismatch in do_em_step");
653 s2->flags = s1->flags;
655 if (s2->natoms != s1->natoms)
657 state_change_natoms(s2, s1->natoms);
658 ems2->f.resize(s2->natoms);
660 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
662 s2->cg_gl.resize(s1->cg_gl.size());
665 copy_mat(s1->box, s2->box);
666 /* Copy free energy state */
667 s2->lambda = s1->lambda;
668 copy_mat(s1->box, s2->box);
673 nthreads = gmx_omp_nthreads_get(emntUpdate);
674 #pragma omp parallel num_threads(nthreads)
676 const rvec* x1 = s1->x.rvec_array();
677 rvec* x2 = s2->x.rvec_array();
678 const rvec* f = as_rvec_array(force.unpaddedArrayRef().data());
681 #pragma omp for schedule(static) nowait
682 for (int i = start; i < end; i++)
690 for (int m = 0; m < DIM; m++)
692 if (ir->opts.nFreeze[gf][m])
698 x2[i][m] = x1[i][m] + a * f[i][m];
702 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
705 if (s2->flags & enumValueToBitMask(StateEntry::Cgp))
707 /* Copy the CG p vector */
708 const rvec* p1 = s1->cg_p.rvec_array();
709 rvec* p2 = s2->cg_p.rvec_array();
710 #pragma omp for schedule(static) nowait
711 for (int i = start; i < end; i++)
713 // Trivial OpenMP block that does not throw
714 copy_rvec(p1[i], p2[i]);
718 if (DOMAINDECOMP(cr))
720 /* OpenMP does not supported unsigned loop variables */
721 #pragma omp for schedule(static) nowait
722 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
724 s2->cg_gl[i] = s1->cg_gl[i];
729 if (DOMAINDECOMP(cr))
731 s2->ddp_count = s1->ddp_count;
732 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
738 validStep = constr->apply(TRUE,
743 s1->x.arrayRefWithPadding(),
744 s2->x.arrayRefWithPadding(),
747 s2->lambda[FreeEnergyPerturbationCouplingType::Bonded],
749 gmx::ArrayRefWithPadding<RVec>(),
752 gmx::ConstraintVariable::Positions);
756 /* This global reduction will affect performance at high
757 * parallelization, but we can not really avoid it.
758 * But usually EM is not run at high parallelization.
760 int reductionBuffer = static_cast<int>(!validStep);
761 gmx_sumi(1, &reductionBuffer, cr);
762 validStep = (reductionBuffer == 0);
765 // We should move this check to the different minimizers
766 if (!validStep && ir->eI != IntegrationAlgorithm::Steep)
769 "The coordinates could not be constrained. Minimizer '%s' can not handle "
770 "constraint failures, use minimizer '%s' before using '%s'.",
771 enumValueToString(ir->eI),
772 enumValueToString(IntegrationAlgorithm::Steep),
773 enumValueToString(ir->eI));
780 //! Prepare EM for using domain decomposition parallellization
781 static void em_dd_partition_system(FILE* fplog,
782 const gmx::MDLogger& mdlog,
785 const gmx_mtop_t& top_global,
786 const t_inputrec* ir,
787 gmx::ImdSession* imdSession,
791 gmx::MDAtoms* mdAtoms,
793 VirtualSitesHandler* vsite,
794 gmx::Constraints* constr,
796 gmx_wallcycle* wcycle)
798 /* Repartition the domain decomposition */
799 dd_partition_system(fplog,
820 dd_store_state(*cr->dd, &ems->s);
826 //! Copy coordinates, OpenMP parallelized, from \p refCoords to coords
827 void setCoordinates(std::vector<RVec>* coords, ArrayRef<const RVec> refCoords)
829 coords->resize(refCoords.size());
831 const int gmx_unused nthreads = gmx_omp_nthreads_get(emntUpdate);
832 #pragma omp parallel for num_threads(nthreads) schedule(static)
833 for (int i = 0; i < ssize(refCoords); i++)
835 (*coords)[i] = refCoords[i];
839 //! Returns the maximum difference an atom moved between two coordinate sets, over all ranks
840 real maxCoordinateDifference(ArrayRef<const RVec> coords1, ArrayRef<const RVec> coords2, MPI_Comm mpiCommMyGroup)
842 GMX_RELEASE_ASSERT(coords1.size() == coords2.size(), "Coordinate counts should match");
844 real maxDiffSquared = 0;
846 const int gmx_unused nthreads = gmx_omp_nthreads_get(emntUpdate);
847 #pragma omp parallel for reduction(max : maxDiffSquared) num_threads(nthreads) schedule(static)
848 for (int i = 0; i < ssize(coords1); i++)
850 maxDiffSquared = std::max(maxDiffSquared, gmx::norm2(coords1[i] - coords2[i]));
855 if (mpiCommMyGroup != MPI_COMM_NULL)
857 MPI_Comm_size(mpiCommMyGroup, &numRanks);
861 real maxDiffSquaredReduced;
863 &maxDiffSquared, &maxDiffSquaredReduced, 1, GMX_DOUBLE ? MPI_DOUBLE : MPI_FLOAT, MPI_MAX, mpiCommMyGroup);
864 maxDiffSquared = maxDiffSquaredReduced;
867 GMX_UNUSED_VALUE(mpiCommMyGroup);
870 return std::sqrt(maxDiffSquared);
873 /*! \brief Class to handle the work of setting and doing an energy evaluation.
875 * This class is a mere aggregate of parameters to pass to evaluate an
876 * energy, so that future changes to names and types of them consume
877 * less time when refactoring other code.
879 * Aggregate initialization is used, for which the chief risk is that
880 * if a member is added at the end and not all initializer lists are
881 * updated, then the member will be value initialized, which will
882 * typically mean initialization to zero.
884 * Use a braced initializer list to construct one of these. */
885 class EnergyEvaluator
888 /*! \brief Evaluates an energy on the state in \c ems.
890 * \todo In practice, the same objects mu_tot, vir, and pres
891 * are always passed to this function, so we would rather have
892 * them as data members. However, their C-array types are
893 * unsuited for aggregate initialization. When the types
894 * improve, the call signature of this method can be reduced.
896 void run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst);
897 //! Handles logging (deprecated).
900 const gmx::MDLogger& mdlog;
901 //! Handles communication.
903 //! Coordinates multi-simulations.
904 const gmx_multisim_t* ms;
905 //! Holds the simulation topology.
906 const gmx_mtop_t& top_global;
907 //! Holds the domain topology.
909 //! User input options.
910 const t_inputrec* inputrec;
911 //! The Interactive Molecular Dynamics session.
912 gmx::ImdSession* imdSession;
913 //! The pull work object.
915 //! Manages flop accounting.
917 //! Manages wall cycle accounting.
918 gmx_wallcycle* wcycle;
919 //! Coordinates global reduction.
920 gmx_global_stat_t gstat;
921 //! Handles virtual sites.
922 VirtualSitesHandler* vsite;
923 //! Handles constraints.
924 gmx::Constraints* constr;
925 //! Per-atom data for this domain.
926 gmx::MDAtoms* mdAtoms;
927 //! Handles how to calculate the forces.
929 //! Schedule of force-calculation work each step for this task.
930 MdrunScheduleWorkload* runScheduleWork;
931 //! Stores the computed energies.
932 gmx_enerdata_t* enerd;
933 //! The DD partitioning count at which the pair list was generated
934 int ddpCountPairSearch;
935 //! The local coordinates that were used for pair searching, stored for computing displacements
936 std::vector<RVec> pairSearchCoordinates;
939 void EnergyEvaluator::run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst)
943 tensor force_vir, shake_vir, ekin;
947 /* Set the time to the initial time, the time does not change during EM */
948 t = inputrec->init_t;
952 vsite->construct(ems->s.x, {}, ems->s.box, gmx::VSiteOperation::Positions);
955 // Compute the buffer size of the pair list
956 const real bufferSize = inputrec->rlist - std::max(inputrec->rcoulomb, inputrec->rvdw);
958 if (bFirst || bufferSize <= 0 || (DOMAINDECOMP(cr) && ems->s.ddp_count != ddpCountPairSearch))
960 /* This is the first state or an old state used before the last ns */
965 // We need to generate a new pairlist when one atom moved more than half the buffer size
966 ArrayRef<const RVec> localCoordinates =
967 ArrayRef<const RVec>(ems->s.x).subArray(0, mdAtoms->mdatoms()->homenr);
968 bNS = 2 * maxCoordinateDifference(pairSearchCoordinates, localCoordinates, cr->mpi_comm_mygroup)
972 if (DOMAINDECOMP(cr) && bNS)
974 /* Repartition the domain decomposition */
975 em_dd_partition_system(
976 fplog, mdlog, count, cr, top_global, inputrec, imdSession, pull_work, ems, top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
977 ddpCountPairSearch = cr->dd->ddp_count;
980 /* Store the local coordinates that will be used in the pair search, after we re-partitioned */
981 if (bufferSize > 0 && bNS)
983 ArrayRef<const RVec> localCoordinates =
984 constArrayRefFromArray(ems->s.x.data(), mdAtoms->mdatoms()->homenr);
985 setCoordinates(&pairSearchCoordinates, localCoordinates);
988 /* Calc force & energy on new trial position */
989 /* do_force always puts the charge groups in the box and shifts again
990 * We do not unshift, so molecules are always whole in congrad.c
1005 ems->s.x.arrayRefWithPadding(),
1018 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY
1019 | (bNS ? GMX_FORCE_NS : 0),
1020 DDBalanceRegionHandler(cr));
1022 /* Clear the unused shake virial and pressure */
1023 clear_mat(shake_vir);
1026 /* Communicate stuff when parallel */
1027 if (PAR(cr) && inputrec->eI != IntegrationAlgorithm::NM)
1029 wallcycle_start(wcycle, WallCycleCounter::MoveE);
1038 gmx::ArrayRef<real>{},
1040 std::vector<real>(1, terminate),
1042 CGLO_ENERGY | CGLO_PRESSURE | CGLO_CONSTRAINT);
1044 wallcycle_stop(wcycle, WallCycleCounter::MoveE);
1047 if (fr->dispersionCorrection)
1049 /* Calculate long range corrections to pressure and energy */
1050 const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
1051 ems->s.box, ems->s.lambda[FreeEnergyPerturbationCouplingType::Vdw]);
1053 enerd->term[F_DISPCORR] = correction.energy;
1054 enerd->term[F_EPOT] += correction.energy;
1055 enerd->term[F_PRES] += correction.pressure;
1056 enerd->term[F_DVDL] += correction.dvdl;
1060 enerd->term[F_DISPCORR] = 0;
1063 ems->epot = enerd->term[F_EPOT];
1067 /* Project out the constraint components of the force */
1068 bool needsLogging = false;
1069 bool computeEnergy = false;
1070 bool computeVirial = true;
1072 auto f = ems->f.view().forceWithPadding();
1073 constr->apply(needsLogging,
1078 ems->s.x.arrayRefWithPadding(),
1080 f.unpaddedArrayRef(),
1082 ems->s.lambda[FreeEnergyPerturbationCouplingType::Bonded],
1084 gmx::ArrayRefWithPadding<RVec>(),
1087 gmx::ConstraintVariable::ForceDispl);
1088 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1089 m_add(force_vir, shake_vir, vir);
1093 copy_mat(force_vir, vir);
1097 enerd->term[F_PRES] = calc_pres(fr->pbcType, inputrec->nwall, ems->s.box, ekin, vir, pres);
1099 if (inputrec->efep != FreeEnergyPerturbationType::No)
1101 accumulateKineticLambdaComponents(enerd, ems->s.lambda, *inputrec->fepvals);
1104 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
1106 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
1112 //! Parallel utility summing energies and forces
1113 static double reorder_partsum(const t_commrec* cr,
1114 const t_grpopts* opts,
1115 const gmx_mtop_t& top_global,
1116 const em_state_t* s_min,
1117 const em_state_t* s_b)
1121 fprintf(debug, "Doing reorder_partsum\n");
1124 auto fm = s_min->f.view().force();
1125 auto fb = s_b->f.view().force();
1127 /* Collect fm in a global vector fmg.
1128 * This conflicts with the spirit of domain decomposition,
1129 * but to fully optimize this a much more complicated algorithm is required.
1131 const int natoms = top_global.natoms;
1135 gmx::ArrayRef<const int> indicesMin = s_min->s.cg_gl;
1137 for (int a : indicesMin)
1139 copy_rvec(fm[i], fmg[a]);
1142 gmx_sum(top_global.natoms * 3, fmg[0], cr);
1144 /* Now we will determine the part of the sum for the cgs in state s_b */
1145 gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
1150 gmx::ArrayRef<const unsigned char> grpnrFREEZE =
1151 top_global.groups.groupNumbers[SimulationAtomGroupType::Freeze];
1152 for (int a : indicesB)
1154 if (!grpnrFREEZE.empty())
1156 gf = grpnrFREEZE[i];
1158 for (int m = 0; m < DIM; m++)
1160 if (!opts->nFreeze[gf][m])
1162 partsum += (fb[i][m] - fmg[a][m]) * fb[i][m];
1173 //! Print some stuff, like beta, whatever that means.
1174 static real pr_beta(const t_commrec* cr,
1175 const t_grpopts* opts,
1177 const gmx_mtop_t& top_global,
1178 const em_state_t* s_min,
1179 const em_state_t* s_b)
1183 /* This is just the classical Polak-Ribiere calculation of beta;
1184 * it looks a bit complicated since we take freeze groups into account,
1185 * and might have to sum it in parallel runs.
1188 if (!DOMAINDECOMP(cr)
1189 || (s_min->s.ddp_count == cr->dd->ddp_count && s_b->s.ddp_count == cr->dd->ddp_count))
1191 auto fm = s_min->f.view().force();
1192 auto fb = s_b->f.view().force();
1195 /* This part of code can be incorrect with DD,
1196 * since the atom ordering in s_b and s_min might differ.
1198 for (int i = 0; i < mdatoms->homenr; i++)
1200 if (mdatoms->cFREEZE)
1202 gf = mdatoms->cFREEZE[i];
1204 for (int m = 0; m < DIM; m++)
1206 if (!opts->nFreeze[gf][m])
1208 sum += (fb[i][m] - fm[i][m]) * fb[i][m];
1215 /* We need to reorder cgs while summing */
1216 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1220 gmx_sumd(1, &sum, cr);
1223 return sum / gmx::square(s_min->fnorm);
1229 void LegacySimulator::do_cg()
1231 const char* CG = "Polak-Ribiere Conjugate Gradients";
1233 gmx_localtop_t top(top_global.ffparams);
1234 gmx_global_stat_t gstat;
1235 double tmp, minstep;
1237 real a, b, c, beta = 0.0;
1240 gmx_bool converged, foundlower;
1241 rvec mu_tot = { 0 };
1242 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1244 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1245 int m, step, nminstep;
1246 auto mdatoms = mdAtoms->mdatoms();
1251 "Note that activating conjugate gradient energy minimization via the "
1252 "integrator .mdp option and the command gmx mdrun may "
1253 "be available in a different form in a future version of GROMACS, "
1254 "e.g. gmx minimize and an .mdp option.");
1260 // In CG, the state is extended with a search direction
1261 state_global->flags |= enumValueToBitMask(StateEntry::Cgp);
1263 // Ensure the extra per-atom state array gets allocated
1264 state_change_natoms(state_global, state_global->natoms);
1266 // Initialize the search direction to zero
1267 for (RVec& cg_p : state_global->cg_p)
1273 /* Create 4 states on the stack and extract pointers that we will swap */
1274 em_state_t s0{}, s1{}, s2{}, s3{};
1275 em_state_t* s_min = &s0;
1276 em_state_t* s_a = &s1;
1277 em_state_t* s_b = &s2;
1278 em_state_t* s_c = &s3;
1280 /* Init em and store the local state in s_min */
1299 const bool simulationsShareState = false;
1300 gmx_mdoutf* outf = init_mdoutf(fplog,
1311 StartingBehavior::NewSimulation,
1312 simulationsShareState,
1314 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
1320 StartingBehavior::NewSimulation,
1321 simulationsShareState,
1322 mdModulesNotifiers);
1324 /* Print to log file */
1325 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1327 /* Max number of steps */
1328 number_steps = inputrec->nsteps;
1332 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1336 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1339 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global,
1340 &top, inputrec, imdSession, pull_work, nrnb,
1341 wcycle, gstat, vsite, constr, mdAtoms,
1342 fr, runScheduleWork, enerd, -1, {} };
1343 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1344 /* do_force always puts the charge groups in the box and shifts again
1345 * We do not unshift, so molecules are always whole in congrad.c
1347 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1351 /* Copy stuff to the energy bin for easy printing etc. */
1352 matrix nullBox = {};
1353 energyOutput.addDataAtEnergyStep(false,
1355 static_cast<double>(step),
1371 EnergyOutput::printHeader(fplog, step, step);
1372 energyOutput.printStepToEnergyFile(
1373 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
1376 /* Estimate/guess the initial stepsize */
1377 stepsize = inputrec->em_stepsize / s_min->fnorm;
1381 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1382 fprintf(stderr, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1383 fprintf(stderr, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1384 fprintf(stderr, "\n");
1385 /* and copy to the log file too... */
1386 fprintf(fplog, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1387 fprintf(fplog, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1388 fprintf(fplog, "\n");
1390 /* Start the loop over CG steps.
1391 * Each successful step is counted, and we continue until
1392 * we either converge or reach the max number of steps.
1395 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1398 /* start taking steps in a new direction
1399 * First time we enter the routine, beta=0, and the direction is
1400 * simply the negative gradient.
1403 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1404 gmx::ArrayRef<gmx::RVec> pm = s_min->s.cg_p;
1405 gmx::ArrayRef<const gmx::RVec> sfm = s_min->f.view().force();
1408 for (int i = 0; i < mdatoms->homenr; i++)
1410 if (mdatoms->cFREEZE)
1412 gf = mdatoms->cFREEZE[i];
1414 for (m = 0; m < DIM; m++)
1416 if (!inputrec->opts.nFreeze[gf][m])
1418 pm[i][m] = sfm[i][m] + beta * pm[i][m];
1419 gpa -= pm[i][m] * sfm[i][m];
1420 /* f is negative gradient, thus the sign */
1429 /* Sum the gradient along the line across CPUs */
1432 gmx_sumd(1, &gpa, cr);
1435 /* Calculate the norm of the search vector */
1436 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1438 /* Just in case stepsize reaches zero due to numerical precision... */
1441 stepsize = inputrec->em_stepsize / pnorm;
1445 * Double check the value of the derivative in the search direction.
1446 * If it is positive it must be due to the old information in the
1447 * CG formula, so just remove that and start over with beta=0.
1448 * This corresponds to a steepest descent step.
1453 step--; /* Don't count this step since we are restarting */
1454 continue; /* Go back to the beginning of the big for-loop */
1457 /* Calculate minimum allowed stepsize, before the average (norm)
1458 * relative change in coordinate is smaller than precision
1461 auto s_min_x = makeArrayRef(s_min->s.x);
1462 for (int i = 0; i < mdatoms->homenr; i++)
1464 for (m = 0; m < DIM; m++)
1466 tmp = fabs(s_min_x[i][m]);
1471 tmp = pm[i][m] / tmp;
1472 minstep += tmp * tmp;
1475 /* Add up from all CPUs */
1478 gmx_sumd(1, &minstep, cr);
1481 minstep = GMX_REAL_EPS / sqrt(minstep / (3 * top_global.natoms));
1483 if (stepsize < minstep)
1489 /* Write coordinates if necessary */
1490 do_x = do_per_step(step, inputrec->nstxout);
1491 do_f = do_per_step(step, inputrec->nstfout);
1494 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, step, s_min, state_global, observablesHistory);
1496 /* Take a step downhill.
1497 * In theory, we should minimize the function along this direction.
1498 * That is quite possible, but it turns out to take 5-10 function evaluations
1499 * for each line. However, we dont really need to find the exact minimum -
1500 * it is much better to start a new CG step in a modified direction as soon
1501 * as we are close to it. This will save a lot of energy evaluations.
1503 * In practice, we just try to take a single step.
1504 * If it worked (i.e. lowered the energy), we increase the stepsize but
1505 * the continue straight to the next CG step without trying to find any minimum.
1506 * If it didn't work (higher energy), there must be a minimum somewhere between
1507 * the old position and the new one.
1509 * Due to the finite numerical accuracy, it turns out that it is a good idea
1510 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1511 * This leads to lower final energies in the tests I've done. / Erik
1513 s_a->epot = s_min->epot;
1515 c = a + stepsize; /* reference position along line is zero */
1517 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1519 em_dd_partition_system(fplog,
1537 /* Take a trial step (new coords in s_c) */
1538 do_em_step(cr, inputrec, mdatoms, s_min, c, s_min->s.cg_p.constArrayRefWithPadding(), s_c, constr, -1);
1541 /* Calculate energy for the trial step */
1542 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1544 /* Calc derivative along line */
1545 const rvec* pc = s_c->s.cg_p.rvec_array();
1546 gmx::ArrayRef<const gmx::RVec> sfc = s_c->f.view().force();
1548 for (int i = 0; i < mdatoms->homenr; i++)
1550 for (m = 0; m < DIM; m++)
1552 gpc -= pc[i][m] * sfc[i][m]; /* f is negative gradient, thus the sign */
1555 /* Sum the gradient along the line across CPUs */
1558 gmx_sumd(1, &gpc, cr);
1561 /* This is the max amount of increase in energy we tolerate */
1562 tmp = std::sqrt(GMX_REAL_EPS) * fabs(s_a->epot);
1564 /* Accept the step if the energy is lower, or if it is not significantly higher
1565 * and the line derivative is still negative.
1567 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1570 /* Great, we found a better energy. Increase step for next iteration
1571 * if we are still going down, decrease it otherwise
1575 stepsize *= 1.618034; /* The golden section */
1579 stepsize *= 0.618034; /* 1/golden section */
1584 /* New energy is the same or higher. We will have to do some work
1585 * to find a smaller value in the interval. Take smaller step next time!
1588 stepsize *= 0.618034;
1592 /* OK, if we didn't find a lower value we will have to locate one now - there must
1593 * be one in the interval [a=0,c].
1594 * The same thing is valid here, though: Don't spend dozens of iterations to find
1595 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1596 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1598 * I also have a safeguard for potentially really pathological functions so we never
1599 * take more than 20 steps before we give up ...
1601 * If we already found a lower value we just skip this step and continue to the update.
1610 /* Select a new trial point.
1611 * If the derivatives at points a & c have different sign we interpolate to zero,
1612 * otherwise just do a bisection.
1614 if (gpa < 0 && gpc > 0)
1616 b = a + gpa * (a - c) / (gpc - gpa);
1623 /* safeguard if interpolation close to machine accuracy causes errors:
1624 * never go outside the interval
1626 if (b <= a || b >= c)
1631 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1633 /* Reload the old state */
1634 em_dd_partition_system(fplog,
1652 /* Take a trial step to this new point - new coords in s_b */
1653 do_em_step(cr, inputrec, mdatoms, s_min, b, s_min->s.cg_p.constArrayRefWithPadding(), s_b, constr, -1);
1656 /* Calculate energy for the trial step */
1657 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1659 /* p does not change within a step, but since the domain decomposition
1660 * might change, we have to use cg_p of s_b here.
1662 const rvec* pb = s_b->s.cg_p.rvec_array();
1663 gmx::ArrayRef<const gmx::RVec> sfb = s_b->f.view().force();
1665 for (int i = 0; i < mdatoms->homenr; i++)
1667 for (m = 0; m < DIM; m++)
1669 gpb -= pb[i][m] * sfb[i][m]; /* f is negative gradient, thus the sign */
1672 /* Sum the gradient along the line across CPUs */
1675 gmx_sumd(1, &gpb, cr);
1680 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a->epot, s_b->epot, s_c->epot, gpb);
1683 epot_repl = s_b->epot;
1685 /* Keep one of the intervals based on the value of the derivative at the new point */
1688 /* Replace c endpoint with b */
1689 swap_em_state(&s_b, &s_c);
1695 /* Replace a endpoint with b */
1696 swap_em_state(&s_b, &s_a);
1702 * Stop search as soon as we find a value smaller than the endpoints.
1703 * Never run more than 20 steps, no matter what.
1706 } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) && (nminstep < 20));
1708 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot) * GMX_REAL_EPS || nminstep >= 20)
1710 /* OK. We couldn't find a significantly lower energy.
1711 * If beta==0 this was steepest descent, and then we give up.
1712 * If not, set beta=0 and restart with steepest descent before quitting.
1722 /* Reset memory before giving up */
1728 /* Select min energy state of A & C, put the best in B.
1730 if (s_c->epot < s_a->epot)
1734 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c->epot, s_a->epot);
1736 swap_em_state(&s_b, &s_c);
1743 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a->epot, s_c->epot);
1745 swap_em_state(&s_b, &s_a);
1753 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n", s_c->epot);
1755 swap_em_state(&s_b, &s_c);
1759 /* new search direction */
1760 /* beta = 0 means forget all memory and restart with steepest descents. */
1761 if (nstcg && ((step % nstcg) == 0))
1767 /* s_min->fnorm cannot be zero, because then we would have converged
1771 /* Polak-Ribiere update.
1772 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1774 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1776 /* Limit beta to prevent oscillations */
1777 if (fabs(beta) > 5.0)
1783 /* update positions */
1784 swap_em_state(&s_min, &s_b);
1787 /* Print it if necessary */
1790 if (mdrunOptions.verbose)
1792 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1794 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1797 s_min->fnorm / sqrtNumAtoms,
1802 /* Store the new (lower) energies */
1803 matrix nullBox = {};
1804 energyOutput.addDataAtEnergyStep(false,
1806 static_cast<double>(step),
1822 do_log = do_per_step(step, inputrec->nstlog);
1823 do_ene = do_per_step(step, inputrec->nstenergy);
1825 imdSession->fillEnergyRecord(step, TRUE);
1829 EnergyOutput::printHeader(fplog, step, step);
1831 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1835 do_log ? fplog : nullptr,
1842 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1843 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x, 0))
1845 imdSession->sendPositionsAndEnergies();
1848 /* Stop when the maximum force lies below tolerance.
1849 * If we have reached machine precision, converged is already set to true.
1851 converged = converged || (s_min->fmax < inputrec->em_tol);
1853 } /* End of the loop */
1857 step--; /* we never took that last step in this case */
1859 if (s_min->fmax > inputrec->em_tol)
1863 warn_step(fplog, inputrec->em_tol, s_min->fmax, step - 1 == number_steps, FALSE);
1870 /* If we printed energy and/or logfile last step (which was the last step)
1871 * we don't have to do it again, but otherwise print the final values.
1875 /* Write final value to log since we didn't do anything the last step */
1876 EnergyOutput::printHeader(fplog, step, step);
1878 if (!do_ene || !do_log)
1880 /* Write final energy file entries */
1881 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1885 !do_log ? fplog : nullptr,
1893 /* Print some stuff... */
1896 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1900 * For accurate normal mode calculation it is imperative that we
1901 * store the last conformation into the full precision binary trajectory.
1903 * However, we should only do it if we did NOT already write this step
1904 * above (which we did if do_x or do_f was true).
1906 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1907 * in the trajectory with the same step number.
1909 do_x = !do_per_step(step, inputrec->nstxout);
1910 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1913 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, s_min, state_global, observablesHistory);
1918 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1919 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1920 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1922 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1925 finish_em(cr, outf, walltime_accounting, wcycle);
1927 /* To print the actual number of steps we needed somewhere */
1928 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1932 void LegacySimulator::do_lbfgs()
1934 static const char* LBFGS = "Low-Memory BFGS Minimizer";
1936 gmx_localtop_t top(top_global.ffparams);
1937 gmx_global_stat_t gstat;
1938 auto mdatoms = mdAtoms->mdatoms();
1943 "Note that activating L-BFGS energy minimization via the "
1944 "integrator .mdp option and the command gmx mdrun may "
1945 "be available in a different form in a future version of GROMACS, "
1946 "e.g. gmx minimize and an .mdp option.");
1950 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1953 if (nullptr != constr)
1957 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1958 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1961 const int n = 3 * state_global->natoms;
1962 const int nmaxcorr = inputrec->nbfgscorr;
1964 std::vector<real> p(n);
1965 std::vector<real> rho(nmaxcorr);
1966 std::vector<real> alpha(nmaxcorr);
1968 std::vector<std::vector<real>> dx(nmaxcorr);
1969 for (auto& dxCorr : dx)
1974 std::vector<std::vector<real>> dg(nmaxcorr);
1975 for (auto& dgCorr : dg)
2002 const bool simulationsShareState = false;
2003 gmx_mdoutf* outf = init_mdoutf(fplog,
2014 StartingBehavior::NewSimulation,
2015 simulationsShareState,
2017 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
2023 StartingBehavior::NewSimulation,
2024 simulationsShareState,
2025 mdModulesNotifiers);
2027 const int start = 0;
2028 const int end = mdatoms->homenr;
2030 /* We need 4 working states */
2031 em_state_t s0{}, s1{}, s2{}, s3{};
2032 em_state_t* sa = &s0;
2033 em_state_t* sb = &s1;
2034 em_state_t* sc = &s2;
2035 em_state_t* last = &s3;
2036 /* Initialize by copying the state from ems (we could skip x and f here) */
2041 /* Print to log file */
2042 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
2044 /* Max number of steps */
2045 const int number_steps = inputrec->nsteps;
2047 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
2048 std::vector<bool> frozen(n);
2050 for (int i = start; i < end; i++)
2052 if (mdatoms->cFREEZE)
2054 gf = mdatoms->cFREEZE[i];
2056 for (int m = 0; m < DIM; m++)
2058 frozen[3 * i + m] = (inputrec->opts.nFreeze[gf][m] != 0);
2063 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
2067 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
2072 vsite->construct(state_global->x, {}, state_global->box, VSiteOperation::Positions);
2075 /* Call the force routine and some auxiliary (neighboursearching etc.) */
2076 /* do_force always puts the charge groups in the box and shifts again
2077 * We do not unshift, so molecules are always whole
2080 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2081 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2082 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2086 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
2090 /* Copy stuff to the energy bin for easy printing etc. */
2091 matrix nullBox = {};
2092 energyOutput.addDataAtEnergyStep(false,
2094 static_cast<double>(step),
2110 EnergyOutput::printHeader(fplog, step, step);
2111 energyOutput.printStepToEnergyFile(
2112 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
2115 /* Set the initial step.
2116 * since it will be multiplied by the non-normalized search direction
2117 * vector (force vector the first time), we scale it by the
2118 * norm of the force.
2123 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2124 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2125 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2126 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2127 fprintf(stderr, "\n");
2128 /* and copy to the log file too... */
2129 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2130 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2131 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2132 fprintf(fplog, "\n");
2135 // Point is an index to the memory of search directions, where 0 is the first one.
2138 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
2139 real* fInit = static_cast<real*>(ems.f.view().force().data()[0]);
2140 for (int i = 0; i < n; i++)
2144 dx[point][i] = fInit[i]; /* Initial search direction */
2152 // Stepsize will be modified during the search, and actually it is not critical
2153 // (the main efficiency in the algorithm comes from changing directions), but
2154 // we still need an initial value, so estimate it as the inverse of the norm
2155 // so we take small steps where the potential fluctuates a lot.
2156 double stepsize = 1.0 / ems.fnorm;
2158 /* Start the loop over BFGS steps.
2159 * Each successful step is counted, and we continue until
2160 * we either converge or reach the max number of steps.
2168 /* Set the gradient from the force */
2169 bool converged = false;
2170 for (int step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
2173 /* Write coordinates if necessary */
2174 const bool do_x = do_per_step(step, inputrec->nstxout);
2175 const bool do_f = do_per_step(step, inputrec->nstfout);
2180 mdof_flags |= MDOF_X;
2185 mdof_flags |= MDOF_F;
2190 mdof_flags |= MDOF_IMD;
2193 gmx::WriteCheckpointDataHolder checkpointDataHolder;
2194 mdoutf_write_to_trajectory_files(fplog,
2200 static_cast<real>(step),
2204 ems.f.view().force(),
2205 &checkpointDataHolder);
2207 /* Do the linesearching in the direction dx[point][0..(n-1)] */
2209 /* make s a pointer to current search direction - point=0 first time we get here */
2210 gmx::ArrayRef<const real> s = dx[point];
2212 const real* xx = static_cast<real*>(ems.s.x.rvec_array()[0]);
2213 const real* ff = static_cast<real*>(ems.f.view().force().data()[0]);
2215 // calculate line gradient in position A
2217 for (int i = 0; i < n; i++)
2219 gpa -= s[i] * ff[i];
2222 /* Calculate minimum allowed stepsize along the line, before the average (norm)
2223 * relative change in coordinate is smaller than precision
2226 for (int i = 0; i < n; i++)
2228 double tmp = fabs(xx[i]);
2234 minstep += tmp * tmp;
2236 minstep = GMX_REAL_EPS / sqrt(minstep / n);
2238 if (stepsize < minstep)
2244 // Before taking any steps along the line, store the old position
2246 real* lastx = static_cast<real*>(last->s.x.data()[0]);
2247 real* lastf = static_cast<real*>(last->f.view().force().data()[0]);
2248 const real Epot0 = ems.epot;
2252 /* Take a step downhill.
2253 * In theory, we should find the actual minimum of the function in this
2254 * direction, somewhere along the line.
2255 * That is quite possible, but it turns out to take 5-10 function evaluations
2256 * for each line. However, we dont really need to find the exact minimum -
2257 * it is much better to start a new BFGS step in a modified direction as soon
2258 * as we are close to it. This will save a lot of energy evaluations.
2260 * In practice, we just try to take a single step.
2261 * If it worked (i.e. lowered the energy), we increase the stepsize but
2262 * continue straight to the next BFGS step without trying to find any minimum,
2263 * i.e. we change the search direction too. If the line was smooth, it is
2264 * likely we are in a smooth region, and then it makes sense to take longer
2265 * steps in the modified search direction too.
2267 * If it didn't work (higher energy), there must be a minimum somewhere between
2268 * the old position and the new one. Then we need to start by finding a lower
2269 * value before we change search direction. Since the energy was apparently
2270 * quite rough, we need to decrease the step size.
2272 * Due to the finite numerical accuracy, it turns out that it is a good idea
2273 * to accept a SMALL increase in energy, if the derivative is still downhill.
2274 * This leads to lower final energies in the tests I've done. / Erik
2277 // State "A" is the first position along the line.
2278 // reference position along line is initially zero
2281 // Check stepsize first. We do not allow displacements
2282 // larger than emstep.
2288 // Pick a new position C by adding stepsize to A.
2291 // Calculate what the largest change in any individual coordinate
2292 // would be (translation along line * gradient along line)
2294 for (int i = 0; i < n; i++)
2296 real delta = c * s[i];
2297 if (delta > maxdelta)
2302 // If any displacement is larger than the stepsize limit, reduce the step
2303 if (maxdelta > inputrec->em_stepsize)
2307 } while (maxdelta > inputrec->em_stepsize);
2309 // Take a trial step and move the coordinate array xc[] to position C
2310 real* xc = static_cast<real*>(sc->s.x.rvec_array()[0]);
2311 for (int i = 0; i < n; i++)
2313 xc[i] = lastx[i] + c * s[i];
2317 // Calculate energy for the trial step in position C
2318 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2320 // Calc line gradient in position C
2321 real* fc = static_cast<real*>(sc->f.view().force()[0]);
2323 for (int i = 0; i < n; i++)
2325 gpc -= s[i] * fc[i]; /* f is negative gradient, thus the sign */
2327 /* Sum the gradient along the line across CPUs */
2330 gmx_sumd(1, &gpc, cr);
2333 // This is the max amount of increase in energy we tolerate.
2334 // By allowing VERY small changes (close to numerical precision) we
2335 // frequently find even better (lower) final energies.
2336 double tmp = std::sqrt(GMX_REAL_EPS) * fabs(sa->epot);
2338 // Accept the step if the energy is lower in the new position C (compared to A),
2339 // or if it is not significantly higher and the line derivative is still negative.
2340 bool foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2341 // If true, great, we found a better energy. We no longer try to alter the
2342 // stepsize, but simply accept this new better position. The we select a new
2343 // search direction instead, which will be much more efficient than continuing
2344 // to take smaller steps along a line. Set fnorm based on the new C position,
2345 // which will be used to update the stepsize to 1/fnorm further down.
2347 // If false, the energy is NOT lower in point C, i.e. it will be the same
2348 // or higher than in point A. In this case it is pointless to move to point C,
2349 // so we will have to do more iterations along the same line to find a smaller
2350 // value in the interval [A=0.0,C].
2351 // Here, A is still 0.0, but that will change when we do a search in the interval
2352 // [0.0,C] below. That search we will do by interpolation or bisection rather
2353 // than with the stepsize, so no need to modify it. For the next search direction
2354 // it will be reset to 1/fnorm anyway.
2359 // OK, if we didn't find a lower value we will have to locate one now - there must
2360 // be one in the interval [a,c].
2361 // The same thing is valid here, though: Don't spend dozens of iterations to find
2362 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2363 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2364 // I also have a safeguard for potentially really pathological functions so we never
2365 // take more than 20 steps before we give up.
2366 // If we already found a lower value we just skip this step and continue to the update.
2371 // Select a new trial point B in the interval [A,C].
2372 // If the derivatives at points a & c have different sign we interpolate to zero,
2373 // otherwise just do a bisection since there might be multiple minima/maxima
2374 // inside the interval.
2376 if (gpa < 0 && gpc > 0)
2378 b = a + gpa * (a - c) / (gpc - gpa);
2385 /* safeguard if interpolation close to machine accuracy causes errors:
2386 * never go outside the interval
2388 if (b <= a || b >= c)
2393 // Take a trial step to point B
2394 real* xb = static_cast<real*>(sb->s.x.rvec_array()[0]);
2395 for (int i = 0; i < n; i++)
2397 xb[i] = lastx[i] + b * s[i];
2401 // Calculate energy for the trial step in point B
2402 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2405 // Calculate gradient in point B
2406 real* fb = static_cast<real*>(sb->f.view().force()[0]);
2408 for (int i = 0; i < n; i++)
2410 gpb -= s[i] * fb[i]; /* f is negative gradient, thus the sign */
2412 /* Sum the gradient along the line across CPUs */
2415 gmx_sumd(1, &gpb, cr);
2418 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2419 // at the new point B, and rename the endpoints of this new interval A and C.
2422 /* Replace c endpoint with b */
2424 /* copy state b to c */
2429 /* Replace a endpoint with b */
2431 /* copy state b to a */
2436 * Stop search as soon as we find a value smaller than the endpoints,
2437 * or if the tolerance is below machine precision.
2438 * Never run more than 20 steps, no matter what.
2441 } while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2443 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2445 /* OK. We couldn't find a significantly lower energy.
2446 * If ncorr==0 this was steepest descent, and then we give up.
2447 * If not, reset memory to restart as steepest descent before quitting.
2459 /* Search in gradient direction */
2460 for (int i = 0; i < n; i++)
2462 dx[point][i] = ff[i];
2464 /* Reset stepsize */
2465 stepsize = 1.0 / fnorm;
2470 /* Select min energy state of A & C, put the best in xx/ff/Epot
2472 if (sc->epot < sa->epot)
2493 /* Update the memory information, and calculate a new
2494 * approximation of the inverse hessian
2497 /* Have new data in Epot, xx, ff */
2498 if (ncorr < nmaxcorr)
2503 for (int i = 0; i < n; i++)
2505 dg[point][i] = lastf[i] - ff[i];
2506 dx[point][i] *= step_taken;
2511 for (int i = 0; i < n; i++)
2513 dgdg += dg[point][i] * dg[point][i];
2514 dgdx += dg[point][i] * dx[point][i];
2517 const real diag = dgdx / dgdg;
2519 rho[point] = 1.0 / dgdx;
2522 if (point >= nmaxcorr)
2528 for (int i = 0; i < n; i++)
2535 /* Recursive update. First go back over the memory points */
2536 for (int k = 0; k < ncorr; k++)
2545 for (int i = 0; i < n; i++)
2547 sq += dx[cp][i] * p[i];
2550 alpha[cp] = rho[cp] * sq;
2552 for (int i = 0; i < n; i++)
2554 p[i] -= alpha[cp] * dg[cp][i];
2558 for (int i = 0; i < n; i++)
2563 /* And then go forward again */
2564 for (int k = 0; k < ncorr; k++)
2567 for (int i = 0; i < n; i++)
2569 yr += p[i] * dg[cp][i];
2572 real beta = rho[cp] * yr;
2573 beta = alpha[cp] - beta;
2575 for (int i = 0; i < n; i++)
2577 p[i] += beta * dx[cp][i];
2587 for (int i = 0; i < n; i++)
2591 dx[point][i] = p[i];
2599 /* Print it if necessary */
2602 if (mdrunOptions.verbose)
2604 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2606 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2609 ems.fnorm / sqrtNumAtoms,
2614 /* Store the new (lower) energies */
2615 matrix nullBox = {};
2616 energyOutput.addDataAtEnergyStep(false,
2618 static_cast<double>(step),
2634 do_log = do_per_step(step, inputrec->nstlog);
2635 do_ene = do_per_step(step, inputrec->nstenergy);
2637 imdSession->fillEnergyRecord(step, TRUE);
2641 EnergyOutput::printHeader(fplog, step, step);
2643 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2647 do_log ? fplog : nullptr,
2654 /* Send x and E to IMD client, if bIMD is TRUE. */
2655 if (imdSession->run(step, TRUE, state_global->box, state_global->x, 0) && MASTER(cr))
2657 imdSession->sendPositionsAndEnergies();
2660 // Reset stepsize in we are doing more iterations
2663 /* Stop when the maximum force lies below tolerance.
2664 * If we have reached machine precision, converged is already set to true.
2666 converged = converged || (ems.fmax < inputrec->em_tol);
2668 } /* End of the loop */
2672 step--; /* we never took that last step in this case */
2674 if (ems.fmax > inputrec->em_tol)
2678 warn_step(fplog, inputrec->em_tol, ems.fmax, step - 1 == number_steps, FALSE);
2683 /* If we printed energy and/or logfile last step (which was the last step)
2684 * we don't have to do it again, but otherwise print the final values.
2686 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2688 EnergyOutput::printHeader(fplog, step, step);
2690 if (!do_ene || !do_log) /* Write final energy file entries */
2692 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2696 !do_log ? fplog : nullptr,
2703 /* Print some stuff... */
2706 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2710 * For accurate normal mode calculation it is imperative that we
2711 * store the last conformation into the full precision binary trajectory.
2713 * However, we should only do it if we did NOT already write this step
2714 * above (which we did if do_x or do_f was true).
2716 const bool do_x = !do_per_step(step, inputrec->nstxout);
2717 const bool do_f = !do_per_step(step, inputrec->nstfout);
2719 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, &ems, state_global, observablesHistory);
2723 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2724 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2725 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2727 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2730 finish_em(cr, outf, walltime_accounting, wcycle);
2732 /* To print the actual number of steps we needed somewhere */
2733 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2736 void LegacySimulator::do_steep()
2738 const char* SD = "Steepest Descents";
2739 gmx_localtop_t top(top_global.ffparams);
2740 gmx_global_stat_t gstat;
2743 gmx_bool bDone, bAbort, do_x, do_f;
2745 rvec mu_tot = { 0 };
2748 int steps_accepted = 0;
2749 auto mdatoms = mdAtoms->mdatoms();
2754 "Note that activating steepest-descent energy minimization via the "
2755 "integrator .mdp option and the command gmx mdrun may "
2756 "be available in a different form in a future version of GROMACS, "
2757 "e.g. gmx minimize and an .mdp option.");
2759 /* Create 2 states on the stack and extract pointers that we will swap */
2760 em_state_t s0{}, s1{};
2761 em_state_t* s_min = &s0;
2762 em_state_t* s_try = &s1;
2764 /* Init em and store the local state in s_try */
2783 const bool simulationsShareState = false;
2784 gmx_mdoutf* outf = init_mdoutf(fplog,
2795 StartingBehavior::NewSimulation,
2796 simulationsShareState,
2798 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
2804 StartingBehavior::NewSimulation,
2805 simulationsShareState,
2806 mdModulesNotifiers);
2808 /* Print to log file */
2809 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2811 /* Set variables for stepsize (in nm). This is the largest
2812 * step that we are going to make in any direction.
2814 ustep = inputrec->em_stepsize;
2817 /* Max number of steps */
2818 nsteps = inputrec->nsteps;
2822 /* Print to the screen */
2823 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2827 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2829 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2830 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2831 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2833 /**** HERE STARTS THE LOOP ****
2834 * count is the counter for the number of steps
2835 * bDone will be TRUE when the minimization has converged
2836 * bAbort will be TRUE when nsteps steps have been performed or when
2837 * the stepsize becomes smaller than is reasonable for machine precision
2842 while (!bDone && !bAbort)
2844 bAbort = (nsteps >= 0) && (count == nsteps);
2846 /* set new coordinates, except for first step */
2847 bool validStep = true;
2850 validStep = do_em_step(
2851 cr, inputrec, mdatoms, s_min, stepsize, s_min->f.view().forceWithPadding(), s_try, constr, count);
2856 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2860 // Signal constraint error during stepping with energy=inf
2861 s_try->epot = std::numeric_limits<real>::infinity();
2866 EnergyOutput::printHeader(fplog, count, count);
2871 s_min->epot = s_try->epot;
2874 /* Print it if necessary */
2877 if (mdrunOptions.verbose)
2880 "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2886 ((count == 0) || (s_try->epot < s_min->epot)) ? '\n' : '\r');
2890 if ((count == 0) || (s_try->epot < s_min->epot))
2892 /* Store the new (lower) energies */
2893 matrix nullBox = {};
2894 energyOutput.addDataAtEnergyStep(false,
2896 static_cast<double>(count),
2912 imdSession->fillEnergyRecord(count, TRUE);
2914 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2915 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2916 energyOutput.printStepToEnergyFile(
2917 mdoutf_get_fp_ene(outf), TRUE, do_dr, do_or, fplog, count, count, fr->fcdata.get(), nullptr);
2922 /* Now if the new energy is smaller than the previous...
2923 * or if this is the first step!
2924 * or if we did random steps!
2927 if ((count == 0) || (s_try->epot < s_min->epot))
2931 /* Test whether the convergence criterion is met... */
2932 bDone = (s_try->fmax < inputrec->em_tol);
2934 /* Copy the arrays for force, positions and energy */
2935 /* The 'Min' array always holds the coords and forces of the minimal
2937 swap_em_state(&s_min, &s_try);
2943 /* Write to trn, if necessary */
2944 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2945 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2947 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, count, s_min, state_global, observablesHistory);
2951 /* If energy is not smaller make the step smaller... */
2954 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2956 /* Reload the old state */
2957 em_dd_partition_system(fplog,
2976 // If the force is very small after finishing minimization,
2977 // we risk dividing by zero when calculating the step size.
2978 // So we check first if the minimization has stopped before
2979 // trying to obtain a new step size.
2982 /* Determine new step */
2983 stepsize = ustep / s_min->fmax;
2986 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2988 if (count == nsteps || ustep < 1e-12)
2990 if (count == nsteps || ustep < 1e-6)
2995 warn_step(fplog, inputrec->em_tol, s_min->fmax, count == nsteps, constr != nullptr);
3000 /* Send IMD energies and positions, if bIMD is TRUE. */
3001 if (imdSession->run(count,
3003 MASTER(cr) ? state_global->box : nullptr,
3004 MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>(),
3008 imdSession->sendPositionsAndEnergies();
3012 } /* End of the loop */
3014 /* Print some data... */
3017 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
3019 write_em_traj(fplog,
3023 inputrec->nstfout != 0,
3024 ftp2fn(efSTO, nfile, fnm),
3030 observablesHistory);
3034 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
3036 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
3037 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
3040 finish_em(cr, outf, walltime_accounting, wcycle);
3042 /* To print the actual number of steps we needed somewhere */
3044 // TODO: Avoid changing inputrec (#3854)
3045 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
3046 nonConstInputrec->nsteps = count;
3049 walltime_accounting_set_nsteps_done(walltime_accounting, count);
3052 void LegacySimulator::do_nm()
3054 const char* NM = "Normal Mode Analysis";
3056 gmx_localtop_t top(top_global.ffparams);
3057 gmx_global_stat_t gstat;
3059 rvec mu_tot = { 0 };
3061 gmx_bool bSparse; /* use sparse matrix storage format */
3063 gmx_sparsematrix_t* sparse_matrix = nullptr;
3064 real* full_matrix = nullptr;
3066 /* added with respect to mdrun */
3068 real der_range = 10.0 * std::sqrt(GMX_REAL_EPS);
3070 bool bIsMaster = MASTER(cr);
3071 auto mdatoms = mdAtoms->mdatoms();
3076 "Note that activating normal-mode analysis via the integrator "
3077 ".mdp option and the command gmx mdrun may "
3078 "be available in a different form in a future version of GROMACS, "
3079 "e.g. gmx normal-modes.");
3081 if (constr != nullptr)
3085 "Constraints present with Normal Mode Analysis, this combination is not supported");
3088 gmx_shellfc_t* shellfc;
3090 em_state_t state_work{};
3092 /* Init em and store the local state in state_minimum */
3111 const bool simulationsShareState = false;
3112 gmx_mdoutf* outf = init_mdoutf(fplog,
3123 StartingBehavior::NewSimulation,
3124 simulationsShareState,
3127 std::vector<int> atom_index = get_atom_index(top_global);
3128 std::vector<gmx::RVec> fneg(atom_index.size(), { 0, 0, 0 });
3129 snew(dfdx, atom_index.size());
3135 "NOTE: This version of GROMACS has been compiled in single precision,\n"
3136 " which MIGHT not be accurate enough for normal mode analysis.\n"
3137 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
3138 " are fairly modest even if you recompile in double precision.\n\n");
3142 /* Check if we can/should use sparse storage format.
3144 * Sparse format is only useful when the Hessian itself is sparse, which it
3145 * will be when we use a cutoff.
3146 * For small systems (n<1000) it is easier to always use full matrix format, though.
3148 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
3150 GMX_LOG(mdlog.warning)
3151 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
3154 else if (atom_index.size() < 1000)
3156 GMX_LOG(mdlog.warning)
3157 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
3163 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
3167 /* Number of dimensions, based on real atoms, that is not vsites or shell */
3168 sz = DIM * atom_index.size();
3170 fprintf(stderr, "Allocating Hessian memory...\n\n");
3174 sparse_matrix = gmx_sparsematrix_init(sz);
3175 sparse_matrix->compressed_symmetric = TRUE;
3179 snew(full_matrix, sz * sz);
3182 /* Write start time and temperature */
3183 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
3185 /* fudge nr of steps to nr of atoms */
3187 // TODO: Avoid changing inputrec (#3854)
3188 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
3189 nonConstInputrec->nsteps = atom_index.size() * 2;
3195 "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
3200 nnodes = cr->nnodes;
3202 /* Make evaluate_energy do a single node force calculation */
3204 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
3205 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
3206 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
3207 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
3208 cr->nnodes = nnodes;
3210 /* if forces are not small, warn user */
3211 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
3213 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
3214 if (state_work.fmax > 1.0e-3)
3216 GMX_LOG(mdlog.warning)
3218 "The force is probably not small enough to "
3219 "ensure that you are at a minimum.\n"
3220 "Be aware that negative eigenvalues may occur\n"
3221 "when the resulting matrix is diagonalized.");
3224 /***********************************************************
3226 * Loop over all pairs in matrix
3228 * do_force called twice. Once with positive and
3229 * once with negative displacement
3231 ************************************************************/
3233 /* Steps are divided one by one over the nodes */
3235 auto state_work_x = makeArrayRef(state_work.s.x);
3236 auto state_work_f = state_work.f.view().force();
3237 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
3239 size_t atom = atom_index[aid];
3240 for (size_t d = 0; d < DIM; d++)
3243 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
3246 x_min = state_work_x[atom][d];
3248 for (unsigned int dx = 0; (dx < 2); dx++)
3252 state_work_x[atom][d] = x_min - der_range;
3256 state_work_x[atom][d] = x_min + der_range;
3259 /* Make evaluate_energy do a single node force calculation */
3263 /* Now is the time to relax the shells */
3264 relax_shell_flexcon(fplog,
3267 mdrunOptions.verbose,
3278 state_work.s.natoms,
3279 state_work.s.x.arrayRefWithPadding(),
3280 state_work.s.v.arrayRefWithPadding(),
3282 state_work.s.lambda,
3284 &state_work.f.view(),
3295 DDBalanceRegionHandler(nullptr));
3301 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid * 2 + dx, FALSE);
3304 cr->nnodes = nnodes;
3308 std::copy(state_work_f.begin(), state_work_f.begin() + atom_index.size(), fneg.begin());
3312 /* x is restored to original */
3313 state_work_x[atom][d] = x_min;
3315 for (size_t j = 0; j < atom_index.size(); j++)
3317 for (size_t k = 0; (k < DIM); k++)
3319 dfdx[j][k] = -(state_work_f[atom_index[j]][k] - fneg[j][k]) / (2 * der_range);
3326 # define mpi_type GMX_MPI_REAL
3327 MPI_Send(dfdx[0], atom_index.size() * DIM, mpi_type, MASTER(cr), cr->nodeid, cr->mpi_comm_mygroup);
3332 for (index node = 0; (node < nnodes && aid + node < ssize(atom_index)); node++)
3338 MPI_Recv(dfdx[0], atom_index.size() * DIM, mpi_type, node, node, cr->mpi_comm_mygroup, &stat);
3343 row = (aid + node) * DIM + d;
3345 for (size_t j = 0; j < atom_index.size(); j++)
3347 for (size_t k = 0; k < DIM; k++)
3353 if (col >= row && dfdx[j][k] != 0.0)
3355 gmx_sparsematrix_increment_value(sparse_matrix, row, col, dfdx[j][k]);
3360 full_matrix[row * sz + col] = dfdx[j][k];
3367 if (mdrunOptions.verbose && fplog)
3372 /* write progress */
3373 if (bIsMaster && mdrunOptions.verbose)
3376 "\rFinished step %d out of %td",
3377 std::min<int>(atom + nnodes, atom_index.size()),
3385 fprintf(stderr, "\n\nWriting Hessian...\n");
3386 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
3389 finish_em(cr, outf, walltime_accounting, wcycle);
3391 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size() * 2);