2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 * \brief This file defines integrators for energy minimization
42 * \author Berk Hess <hess@kth.se>
43 * \author Erik Lindahl <erik@kth.se>
44 * \ingroup module_mdrun
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/collect.h"
59 #include "gromacs/domdec/dlbtiming.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/mdsetup.h"
63 #include "gromacs/domdec/partition.h"
64 #include "gromacs/ewald/pme_pp.h"
65 #include "gromacs/fileio/confio.h"
66 #include "gromacs/fileio/mtxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/imd/imd.h"
70 #include "gromacs/linearalgebra/sparsematrix.h"
71 #include "gromacs/listed_forces/listed_forces.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/mdlib/constr.h"
75 #include "gromacs/mdlib/coupling.h"
76 #include "gromacs/mdlib/dispersioncorrection.h"
77 #include "gromacs/mdlib/ebin.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/energyoutput.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/gmx_omp_nthreads.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/stat.h"
87 #include "gromacs/mdlib/tgroup.h"
88 #include "gromacs/mdlib/trajectory_writing.h"
89 #include "gromacs/mdlib/update.h"
90 #include "gromacs/mdlib/vsite.h"
91 #include "gromacs/mdrunutility/handlerestart.h"
92 #include "gromacs/mdrunutility/printtime.h"
93 #include "gromacs/mdtypes/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_t wcycle,
146 walltime_accounting_start_time(walltime_accounting);
147 wallcycle_start(wcycle, ewcRUN);
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_t wcycle)
154 wallcycle_stop(wcycle, ewcRUN);
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",
240 alg, gmx_step_str(count, buf), ftol);
244 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n", alg, ftol,
245 gmx_step_str(count, buf));
249 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
250 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
251 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm / sqrtNumAtoms);
253 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
254 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
255 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm / sqrtNumAtoms);
259 //! Compute the norm and max of the force array in parallel
260 static void get_f_norm_max(const t_commrec* cr,
263 gmx::ArrayRef<const gmx::RVec> f,
270 int la_max, a_max, start, end, i, m, gf;
272 /* This routine finds the largest force and returns it.
273 * On parallel machines the global max is taken.
279 end = mdatoms->homenr;
280 if (mdatoms->cFREEZE)
282 for (i = start; i < end; i++)
284 gf = mdatoms->cFREEZE[i];
286 for (m = 0; m < DIM; m++)
288 if (!opts->nFreeze[gf][m])
290 fam += gmx::square(f[i][m]);
303 for (i = start; i < end; i++)
315 if (la_max >= 0 && DOMAINDECOMP(cr))
317 a_max = cr->dd->globalAtomIndices[la_max];
325 snew(sum, 2 * cr->nnodes + 1);
326 sum[2 * cr->nodeid] = fmax2;
327 sum[2 * cr->nodeid + 1] = a_max;
328 sum[2 * cr->nnodes] = fnorm2;
329 gmx_sumd(2 * cr->nnodes + 1, sum, cr);
330 fnorm2 = sum[2 * cr->nnodes];
331 /* Determine the global maximum */
332 for (i = 0; i < cr->nnodes; i++)
334 if (sum[2 * i] > fmax2)
337 a_max = gmx::roundToInt(sum[2 * i + 1]);
345 *fnorm = sqrt(fnorm2);
357 //! Compute the norm of the force
358 static void get_state_f_norm_max(const t_commrec* cr, t_grpopts* opts, t_mdatoms* mdatoms, em_state_t* ems)
360 get_f_norm_max(cr, opts, mdatoms, ems->f.view().force(), &ems->fnorm, &ems->fmax, &ems->a_fmax);
363 //! Initialize the energy minimization
364 static void init_em(FILE* fplog,
365 const gmx::MDLogger& mdlog,
369 gmx::ImdSession* imdSession,
371 t_state* state_global,
372 const gmx_mtop_t* top_global,
377 gmx::MDAtoms* mdAtoms,
378 gmx_global_stat_t* gstat,
379 VirtualSitesHandler* vsite,
380 gmx::Constraints* constr,
381 gmx_shellfc_t** shellfc)
387 fprintf(fplog, "Initiating %s\n", title);
392 state_global->ngtc = 0;
394 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
395 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
396 initialize_lambdas(fplog, *ir, MASTER(cr), fep_state, lambda, nullptr);
400 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
402 *shellfc = init_shell_flexcon(stdout, top_global, constr ? constr->numFlexibleConstraints() : 0,
403 ir->nstcalcenergy, DOMAINDECOMP(cr));
407 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI),
408 "This else currently only handles energy minimizers, consider if your algorithm "
409 "needs shell/flexible-constraint support");
411 /* With energy minimization, shells and flexible constraints are
412 * automatically minimized when treated like normal DOFS.
414 if (shellfc != nullptr)
420 if (DOMAINDECOMP(cr))
422 dd_init_local_state(cr->dd, state_global, &ems->s);
424 /* Distribute the charge groups over the nodes from the master node */
425 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
426 imdSession, pull_work, &ems->s, &ems->f, mdAtoms, top, fr, vsite,
427 constr, nrnb, nullptr, FALSE);
428 dd_store_state(cr->dd, &ems->s);
432 state_change_natoms(state_global, state_global->natoms);
433 /* Just copy the state */
434 ems->s = *state_global;
435 state_change_natoms(&ems->s, ems->s.natoms);
437 mdAlgorithmsSetupAtomData(cr, ir, *top_global, top, fr, &ems->f, mdAtoms, constr, vsite,
438 shellfc ? *shellfc : nullptr);
441 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
445 // TODO how should this cross-module support dependency be managed?
446 if (ir->eConstrAlg == econtSHAKE && gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
448 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
449 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
452 if (!ir->bContinuation)
454 /* Constrain the starting coordinates */
455 bool needsLogging = true;
456 bool computeEnergy = true;
457 bool computeVirial = false;
459 constr->apply(needsLogging, computeEnergy, -1, 0, 1.0, ems->s.x.arrayRefWithPadding(),
460 ems->s.x.arrayRefWithPadding(), ArrayRef<RVec>(), ems->s.box,
461 ems->s.lambda[efptFEP], &dvdl_constr, gmx::ArrayRefWithPadding<RVec>(),
462 computeVirial, nullptr, gmx::ConstraintVariable::Positions);
468 *gstat = global_stat_init(ir);
475 calc_shifts(ems->s.box, fr->shift_vec);
478 //! Finalize the minimization
479 static void finish_em(const t_commrec* cr,
481 gmx_walltime_accounting_t walltime_accounting,
482 gmx_wallcycle_t wcycle)
484 if (!thisRankHasDuty(cr, DUTY_PME))
486 /* Tell the PME only node to finish */
487 gmx_pme_send_finish(cr);
492 em_time_end(walltime_accounting, wcycle);
495 //! Swap two different EM states during minimization
496 static void swap_em_state(em_state_t** ems1, em_state_t** ems2)
505 //! Save the EM trajectory
506 static void write_em_traj(FILE* fplog,
512 const gmx_mtop_t* top_global,
516 t_state* state_global,
517 ObservablesHistory* observablesHistory)
523 mdof_flags |= MDOF_X;
527 mdof_flags |= MDOF_F;
530 /* If we want IMD output, set appropriate MDOF flag */
533 mdof_flags |= MDOF_IMD;
536 gmx::WriteCheckpointDataHolder checkpointDataHolder;
537 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global->natoms, step,
538 static_cast<double>(step), &state->s, state_global,
539 observablesHistory, state->f.view().force(), &checkpointDataHolder);
541 if (confout != nullptr)
543 if (DOMAINDECOMP(cr))
545 /* If bX=true, x was collected to state_global in the call above */
548 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
549 dd_collect_vec(cr->dd, state->s.ddp_count, state->s.ddp_count_cg_gl, state->s.cg_gl,
550 state->s.x, globalXRef);
555 /* Copy the local state pointer */
556 state_global = &state->s;
561 if (ir->pbcType != PbcType::No && !ir->bPeriodicMols && DOMAINDECOMP(cr))
563 /* Make molecules whole only for confout writing */
564 do_pbc_mtop(ir->pbcType, state->s.box, top_global, state_global->x.rvec_array());
567 write_sto_conf_mtop(confout, *top_global->name, top_global,
568 state_global->x.rvec_array(), nullptr, ir->pbcType, state->s.box);
573 //! \brief Do one minimization step
575 // \returns true when the step succeeded, false when a constraint error occurred
576 static bool do_em_step(const t_commrec* cr,
581 gmx::ArrayRefWithPadding<const gmx::RVec> force,
583 gmx::Constraints* constr,
590 int nthreads gmx_unused;
592 bool validStep = true;
597 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
599 gmx_incons("state mismatch in do_em_step");
602 s2->flags = s1->flags;
604 if (s2->natoms != s1->natoms)
606 state_change_natoms(s2, s1->natoms);
607 ems2->f.resize(s2->natoms);
609 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
611 s2->cg_gl.resize(s1->cg_gl.size());
614 copy_mat(s1->box, s2->box);
615 /* Copy free energy state */
616 s2->lambda = s1->lambda;
617 copy_mat(s1->box, s2->box);
622 nthreads = gmx_omp_nthreads_get(emntUpdate);
623 #pragma omp parallel num_threads(nthreads)
625 const rvec* x1 = s1->x.rvec_array();
626 rvec* x2 = s2->x.rvec_array();
627 const rvec* f = as_rvec_array(force.unpaddedArrayRef().data());
630 #pragma omp for schedule(static) nowait
631 for (int i = start; i < end; i++)
639 for (int m = 0; m < DIM; m++)
641 if (ir->opts.nFreeze[gf][m])
647 x2[i][m] = x1[i][m] + a * f[i][m];
651 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
654 if (s2->flags & (1 << estCGP))
656 /* Copy the CG p vector */
657 const rvec* p1 = s1->cg_p.rvec_array();
658 rvec* p2 = s2->cg_p.rvec_array();
659 #pragma omp for schedule(static) nowait
660 for (int i = start; i < end; i++)
662 // Trivial OpenMP block that does not throw
663 copy_rvec(p1[i], p2[i]);
667 if (DOMAINDECOMP(cr))
669 /* OpenMP does not supported unsigned loop variables */
670 #pragma omp for schedule(static) nowait
671 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
673 s2->cg_gl[i] = s1->cg_gl[i];
678 if (DOMAINDECOMP(cr))
680 s2->ddp_count = s1->ddp_count;
681 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
687 validStep = constr->apply(
688 TRUE, TRUE, count, 0, 1.0, s1->x.arrayRefWithPadding(), s2->x.arrayRefWithPadding(),
689 ArrayRef<RVec>(), s2->box, s2->lambda[efptBONDED], &dvdl_constr,
690 gmx::ArrayRefWithPadding<RVec>(), false, nullptr, gmx::ConstraintVariable::Positions);
694 /* This global reduction will affect performance at high
695 * parallelization, but we can not really avoid it.
696 * But usually EM is not run at high parallelization.
698 int reductionBuffer = static_cast<int>(!validStep);
699 gmx_sumi(1, &reductionBuffer, cr);
700 validStep = (reductionBuffer == 0);
703 // We should move this check to the different minimizers
704 if (!validStep && ir->eI != eiSteep)
707 "The coordinates could not be constrained. Minimizer '%s' can not handle "
708 "constraint failures, use minimizer '%s' before using '%s'.",
709 EI(ir->eI), EI(eiSteep), EI(ir->eI));
716 //! Prepare EM for using domain decomposition parallellization
717 static void em_dd_partition_system(FILE* fplog,
718 const gmx::MDLogger& mdlog,
721 const gmx_mtop_t* top_global,
723 gmx::ImdSession* imdSession,
727 gmx::MDAtoms* mdAtoms,
729 VirtualSitesHandler* vsite,
730 gmx::Constraints* constr,
732 gmx_wallcycle_t wcycle)
734 /* Repartition the domain decomposition */
735 dd_partition_system(fplog, mdlog, step, cr, FALSE, 1, nullptr, *top_global, ir, imdSession, pull_work,
736 &ems->s, &ems->f, mdAtoms, top, fr, vsite, constr, nrnb, wcycle, FALSE);
737 dd_store_state(cr->dd, &ems->s);
743 /*! \brief Class to handle the work of setting and doing an energy evaluation.
745 * This class is a mere aggregate of parameters to pass to evaluate an
746 * energy, so that future changes to names and types of them consume
747 * less time when refactoring other code.
749 * Aggregate initialization is used, for which the chief risk is that
750 * if a member is added at the end and not all initializer lists are
751 * updated, then the member will be value initialized, which will
752 * typically mean initialization to zero.
754 * Use a braced initializer list to construct one of these. */
755 class EnergyEvaluator
758 /*! \brief Evaluates an energy on the state in \c ems.
760 * \todo In practice, the same objects mu_tot, vir, and pres
761 * are always passed to this function, so we would rather have
762 * them as data members. However, their C-array types are
763 * unsuited for aggregate initialization. When the types
764 * improve, the call signature of this method can be reduced.
766 void run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst);
767 //! Handles logging (deprecated).
770 const gmx::MDLogger& mdlog;
771 //! Handles communication.
773 //! Coordinates multi-simulations.
774 const gmx_multisim_t* ms;
775 //! Holds the simulation topology.
776 const gmx_mtop_t* top_global;
777 //! Holds the domain topology.
779 //! User input options.
780 t_inputrec* inputrec;
781 //! The Interactive Molecular Dynamics session.
782 gmx::ImdSession* imdSession;
783 //! The pull work object.
785 //! Manages flop accounting.
787 //! Manages wall cycle accounting.
788 gmx_wallcycle_t wcycle;
789 //! Coordinates global reduction.
790 gmx_global_stat_t gstat;
791 //! Handles virtual sites.
792 VirtualSitesHandler* vsite;
793 //! Handles constraints.
794 gmx::Constraints* constr;
795 //! Per-atom data for this domain.
796 gmx::MDAtoms* mdAtoms;
797 //! Handles how to calculate the forces.
799 //! Schedule of force-calculation work each step for this task.
800 MdrunScheduleWorkload* runScheduleWork;
801 //! Stores the computed energies.
802 gmx_enerdata_t* enerd;
805 void EnergyEvaluator::run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst)
809 tensor force_vir, shake_vir, ekin;
813 /* Set the time to the initial time, the time does not change during EM */
814 t = inputrec->init_t;
816 if (bFirst || (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
818 /* This is the first state or an old state used before the last ns */
824 if (inputrec->nstlist > 0)
832 vsite->construct(ems->s.x, 1, {}, ems->s.box);
835 if (DOMAINDECOMP(cr) && bNS)
837 /* Repartition the domain decomposition */
838 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession, pull_work,
839 ems, top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
842 /* Calc force & energy on new trial position */
843 /* do_force always puts the charge groups in the box and shifts again
844 * We do not unshift, so molecules are always whole in congrad.c
846 do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession, pull_work, count, nrnb, wcycle,
847 top, ems->s.box, ems->s.x.arrayRefWithPadding(), &ems->s.hist, &ems->f.view(), force_vir,
848 mdAtoms->mdatoms(), enerd, ems->s.lambda, fr, runScheduleWork, vsite, mu_tot, t, nullptr,
849 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY
850 | (bNS ? GMX_FORCE_NS : 0),
851 DDBalanceRegionHandler(cr));
853 /* Clear the unused shake virial and pressure */
854 clear_mat(shake_vir);
857 /* Communicate stuff when parallel */
858 if (PAR(cr) && inputrec->eI != eiNM)
860 wallcycle_start(wcycle, ewcMoveE);
862 global_stat(gstat, cr, enerd, force_vir, shake_vir, inputrec, nullptr, nullptr, nullptr, 1,
863 &terminate, nullptr, FALSE, CGLO_ENERGY | CGLO_PRESSURE | CGLO_CONSTRAINT);
865 wallcycle_stop(wcycle, ewcMoveE);
868 if (fr->dispersionCorrection)
870 /* Calculate long range corrections to pressure and energy */
871 const DispersionCorrection::Correction correction =
872 fr->dispersionCorrection->calculate(ems->s.box, ems->s.lambda[efptVDW]);
874 enerd->term[F_DISPCORR] = correction.energy;
875 enerd->term[F_EPOT] += correction.energy;
876 enerd->term[F_PRES] += correction.pressure;
877 enerd->term[F_DVDL] += correction.dvdl;
881 enerd->term[F_DISPCORR] = 0;
884 ems->epot = enerd->term[F_EPOT];
888 /* Project out the constraint components of the force */
889 bool needsLogging = false;
890 bool computeEnergy = false;
891 bool computeVirial = true;
893 auto f = ems->f.view().forceWithPadding();
894 constr->apply(needsLogging, computeEnergy, count, 0, 1.0, ems->s.x.arrayRefWithPadding(), f,
895 f.unpaddedArrayRef(), ems->s.box, ems->s.lambda[efptBONDED], &dvdl_constr,
896 gmx::ArrayRefWithPadding<RVec>(), computeVirial, shake_vir,
897 gmx::ConstraintVariable::ForceDispl);
898 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
899 m_add(force_vir, shake_vir, vir);
903 copy_mat(force_vir, vir);
907 enerd->term[F_PRES] = calc_pres(fr->pbcType, inputrec->nwall, ems->s.box, ekin, vir, pres);
909 if (inputrec->efep != efepNO)
911 accumulateKineticLambdaComponents(enerd, ems->s.lambda, *inputrec->fepvals);
914 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
916 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
922 //! Parallel utility summing energies and forces
923 static double reorder_partsum(const t_commrec* cr,
925 const gmx_mtop_t* top_global,
926 const em_state_t* s_min,
927 const em_state_t* s_b)
931 fprintf(debug, "Doing reorder_partsum\n");
934 auto fm = s_min->f.view().force();
935 auto fb = s_b->f.view().force();
937 /* Collect fm in a global vector fmg.
938 * This conflicts with the spirit of domain decomposition,
939 * but to fully optimize this a much more complicated algorithm is required.
941 const int natoms = top_global->natoms;
945 gmx::ArrayRef<const int> indicesMin = s_min->s.cg_gl;
947 for (int a : indicesMin)
949 copy_rvec(fm[i], fmg[a]);
952 gmx_sum(top_global->natoms * 3, fmg[0], cr);
954 /* Now we will determine the part of the sum for the cgs in state s_b */
955 gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
960 gmx::ArrayRef<const unsigned char> grpnrFREEZE =
961 top_global->groups.groupNumbers[SimulationAtomGroupType::Freeze];
962 for (int a : indicesB)
964 if (!grpnrFREEZE.empty())
968 for (int m = 0; m < DIM; m++)
970 if (!opts->nFreeze[gf][m])
972 partsum += (fb[i][m] - fmg[a][m]) * fb[i][m];
983 //! Print some stuff, like beta, whatever that means.
984 static real pr_beta(const t_commrec* cr,
987 const gmx_mtop_t* top_global,
988 const em_state_t* s_min,
989 const em_state_t* s_b)
993 /* This is just the classical Polak-Ribiere calculation of beta;
994 * it looks a bit complicated since we take freeze groups into account,
995 * and might have to sum it in parallel runs.
998 if (!DOMAINDECOMP(cr)
999 || (s_min->s.ddp_count == cr->dd->ddp_count && s_b->s.ddp_count == cr->dd->ddp_count))
1001 auto fm = s_min->f.view().force();
1002 auto fb = s_b->f.view().force();
1005 /* This part of code can be incorrect with DD,
1006 * since the atom ordering in s_b and s_min might differ.
1008 for (int i = 0; i < mdatoms->homenr; i++)
1010 if (mdatoms->cFREEZE)
1012 gf = mdatoms->cFREEZE[i];
1014 for (int m = 0; m < DIM; m++)
1016 if (!opts->nFreeze[gf][m])
1018 sum += (fb[i][m] - fm[i][m]) * fb[i][m];
1025 /* We need to reorder cgs while summing */
1026 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1030 gmx_sumd(1, &sum, cr);
1033 return sum / gmx::square(s_min->fnorm);
1039 void LegacySimulator::do_cg()
1041 const char* CG = "Polak-Ribiere Conjugate Gradients";
1043 gmx_localtop_t top(top_global->ffparams);
1044 gmx_global_stat_t gstat;
1045 double tmp, minstep;
1047 real a, b, c, beta = 0.0;
1050 gmx_bool converged, foundlower;
1051 rvec mu_tot = { 0 };
1052 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1054 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1055 int m, step, nminstep;
1056 auto mdatoms = mdAtoms->mdatoms();
1061 "Note that activating conjugate gradient energy minimization via the "
1062 "integrator .mdp option and the command gmx mdrun may "
1063 "be available in a different form in a future version of GROMACS, "
1064 "e.g. gmx minimize and an .mdp option.");
1070 // In CG, the state is extended with a search direction
1071 state_global->flags |= (1 << estCGP);
1073 // Ensure the extra per-atom state array gets allocated
1074 state_change_natoms(state_global, state_global->natoms);
1076 // Initialize the search direction to zero
1077 for (RVec& cg_p : state_global->cg_p)
1083 /* Create 4 states on the stack and extract pointers that we will swap */
1084 em_state_t s0{}, s1{}, s2{}, s3{};
1085 em_state_t* s_min = &s0;
1086 em_state_t* s_a = &s1;
1087 em_state_t* s_b = &s2;
1088 em_state_t* s_c = &s3;
1090 /* Init em and store the local state in s_min */
1091 init_em(fplog, mdlog, CG, cr, inputrec, imdSession, pull_work, state_global, top_global, s_min,
1092 &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, nullptr);
1093 const bool simulationsShareState = false;
1094 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
1095 mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1096 StartingBehavior::NewSimulation, simulationsShareState, ms);
1097 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
1098 false, StartingBehavior::NewSimulation, mdModulesNotifier);
1100 /* Print to log file */
1101 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1103 /* Max number of steps */
1104 number_steps = inputrec->nsteps;
1108 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1112 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1115 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
1116 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1117 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
1118 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1119 /* do_force always puts the charge groups in the box and shifts again
1120 * We do not unshift, so molecules are always whole in congrad.c
1122 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1126 /* Copy stuff to the energy bin for easy printing etc. */
1127 matrix nullBox = {};
1128 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1129 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1130 nullptr, vir, pres, nullptr, mu_tot, constr);
1132 EnergyOutput::printHeader(fplog, step, step);
1133 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step,
1134 step, fr->fcdata.get(), nullptr);
1137 /* Estimate/guess the initial stepsize */
1138 stepsize = inputrec->em_stepsize / s_min->fnorm;
1142 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1143 fprintf(stderr, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1144 fprintf(stderr, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1145 fprintf(stderr, "\n");
1146 /* and copy to the log file too... */
1147 fprintf(fplog, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1148 fprintf(fplog, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1149 fprintf(fplog, "\n");
1151 /* Start the loop over CG steps.
1152 * Each successful step is counted, and we continue until
1153 * we either converge or reach the max number of steps.
1156 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1159 /* start taking steps in a new direction
1160 * First time we enter the routine, beta=0, and the direction is
1161 * simply the negative gradient.
1164 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1165 gmx::ArrayRef<gmx::RVec> pm = s_min->s.cg_p;
1166 gmx::ArrayRef<const gmx::RVec> sfm = s_min->f.view().force();
1169 for (int i = 0; i < mdatoms->homenr; i++)
1171 if (mdatoms->cFREEZE)
1173 gf = mdatoms->cFREEZE[i];
1175 for (m = 0; m < DIM; m++)
1177 if (!inputrec->opts.nFreeze[gf][m])
1179 pm[i][m] = sfm[i][m] + beta * pm[i][m];
1180 gpa -= pm[i][m] * sfm[i][m];
1181 /* f is negative gradient, thus the sign */
1190 /* Sum the gradient along the line across CPUs */
1193 gmx_sumd(1, &gpa, cr);
1196 /* Calculate the norm of the search vector */
1197 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1199 /* Just in case stepsize reaches zero due to numerical precision... */
1202 stepsize = inputrec->em_stepsize / pnorm;
1206 * Double check the value of the derivative in the search direction.
1207 * If it is positive it must be due to the old information in the
1208 * CG formula, so just remove that and start over with beta=0.
1209 * This corresponds to a steepest descent step.
1214 step--; /* Don't count this step since we are restarting */
1215 continue; /* Go back to the beginning of the big for-loop */
1218 /* Calculate minimum allowed stepsize, before the average (norm)
1219 * relative change in coordinate is smaller than precision
1222 auto s_min_x = makeArrayRef(s_min->s.x);
1223 for (int i = 0; i < mdatoms->homenr; i++)
1225 for (m = 0; m < DIM; m++)
1227 tmp = fabs(s_min_x[i][m]);
1232 tmp = pm[i][m] / tmp;
1233 minstep += tmp * tmp;
1236 /* Add up from all CPUs */
1239 gmx_sumd(1, &minstep, cr);
1242 minstep = GMX_REAL_EPS / sqrt(minstep / (3 * top_global->natoms));
1244 if (stepsize < minstep)
1250 /* Write coordinates if necessary */
1251 do_x = do_per_step(step, inputrec->nstxout);
1252 do_f = do_per_step(step, inputrec->nstfout);
1254 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, step, s_min,
1255 state_global, observablesHistory);
1257 /* Take a step downhill.
1258 * In theory, we should minimize the function along this direction.
1259 * That is quite possible, but it turns out to take 5-10 function evaluations
1260 * for each line. However, we dont really need to find the exact minimum -
1261 * it is much better to start a new CG step in a modified direction as soon
1262 * as we are close to it. This will save a lot of energy evaluations.
1264 * In practice, we just try to take a single step.
1265 * If it worked (i.e. lowered the energy), we increase the stepsize but
1266 * the continue straight to the next CG step without trying to find any minimum.
1267 * If it didn't work (higher energy), there must be a minimum somewhere between
1268 * the old position and the new one.
1270 * Due to the finite numerical accuracy, it turns out that it is a good idea
1271 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1272 * This leads to lower final energies in the tests I've done. / Erik
1274 s_a->epot = s_min->epot;
1276 c = a + stepsize; /* reference position along line is zero */
1278 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1280 em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec, imdSession,
1281 pull_work, s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
1284 /* Take a trial step (new coords in s_c) */
1285 do_em_step(cr, inputrec, mdatoms, s_min, c, s_min->s.cg_p.constArrayRefWithPadding(), s_c,
1289 /* Calculate energy for the trial step */
1290 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1292 /* Calc derivative along line */
1293 const rvec* pc = s_c->s.cg_p.rvec_array();
1294 gmx::ArrayRef<const gmx::RVec> sfc = s_c->f.view().force();
1296 for (int i = 0; i < mdatoms->homenr; i++)
1298 for (m = 0; m < DIM; m++)
1300 gpc -= pc[i][m] * sfc[i][m]; /* f is negative gradient, thus the sign */
1303 /* Sum the gradient along the line across CPUs */
1306 gmx_sumd(1, &gpc, cr);
1309 /* This is the max amount of increase in energy we tolerate */
1310 tmp = std::sqrt(GMX_REAL_EPS) * fabs(s_a->epot);
1312 /* Accept the step if the energy is lower, or if it is not significantly higher
1313 * and the line derivative is still negative.
1315 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1318 /* Great, we found a better energy. Increase step for next iteration
1319 * if we are still going down, decrease it otherwise
1323 stepsize *= 1.618034; /* The golden section */
1327 stepsize *= 0.618034; /* 1/golden section */
1332 /* New energy is the same or higher. We will have to do some work
1333 * to find a smaller value in the interval. Take smaller step next time!
1336 stepsize *= 0.618034;
1340 /* OK, if we didn't find a lower value we will have to locate one now - there must
1341 * be one in the interval [a=0,c].
1342 * The same thing is valid here, though: Don't spend dozens of iterations to find
1343 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1344 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1346 * I also have a safeguard for potentially really pathological functions so we never
1347 * take more than 20 steps before we give up ...
1349 * If we already found a lower value we just skip this step and continue to the update.
1358 /* Select a new trial point.
1359 * If the derivatives at points a & c have different sign we interpolate to zero,
1360 * otherwise just do a bisection.
1362 if (gpa < 0 && gpc > 0)
1364 b = a + gpa * (a - c) / (gpc - gpa);
1371 /* safeguard if interpolation close to machine accuracy causes errors:
1372 * never go outside the interval
1374 if (b <= a || b >= c)
1379 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1381 /* Reload the old state */
1382 em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec, imdSession, pull_work,
1383 s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
1386 /* Take a trial step to this new point - new coords in s_b */
1387 do_em_step(cr, inputrec, mdatoms, s_min, b,
1388 s_min->s.cg_p.constArrayRefWithPadding(), s_b, constr, -1);
1391 /* Calculate energy for the trial step */
1392 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1394 /* p does not change within a step, but since the domain decomposition
1395 * might change, we have to use cg_p of s_b here.
1397 const rvec* pb = s_b->s.cg_p.rvec_array();
1398 gmx::ArrayRef<const gmx::RVec> sfb = s_b->f.view().force();
1400 for (int i = 0; i < mdatoms->homenr; i++)
1402 for (m = 0; m < DIM; m++)
1404 gpb -= pb[i][m] * sfb[i][m]; /* f is negative gradient, thus the sign */
1407 /* Sum the gradient along the line across CPUs */
1410 gmx_sumd(1, &gpb, cr);
1415 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a->epot, s_b->epot,
1419 epot_repl = s_b->epot;
1421 /* Keep one of the intervals based on the value of the derivative at the new point */
1424 /* Replace c endpoint with b */
1425 swap_em_state(&s_b, &s_c);
1431 /* Replace a endpoint with b */
1432 swap_em_state(&s_b, &s_a);
1438 * Stop search as soon as we find a value smaller than the endpoints.
1439 * Never run more than 20 steps, no matter what.
1442 } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) && (nminstep < 20));
1444 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot) * GMX_REAL_EPS || nminstep >= 20)
1446 /* OK. We couldn't find a significantly lower energy.
1447 * If beta==0 this was steepest descent, and then we give up.
1448 * If not, set beta=0 and restart with steepest descent before quitting.
1458 /* Reset memory before giving up */
1464 /* Select min energy state of A & C, put the best in B.
1466 if (s_c->epot < s_a->epot)
1470 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c->epot,
1473 swap_em_state(&s_b, &s_c);
1480 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a->epot,
1483 swap_em_state(&s_b, &s_a);
1491 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n", s_c->epot);
1493 swap_em_state(&s_b, &s_c);
1497 /* new search direction */
1498 /* beta = 0 means forget all memory and restart with steepest descents. */
1499 if (nstcg && ((step % nstcg) == 0))
1505 /* s_min->fnorm cannot be zero, because then we would have converged
1509 /* Polak-Ribiere update.
1510 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1512 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1514 /* Limit beta to prevent oscillations */
1515 if (fabs(beta) > 5.0)
1521 /* update positions */
1522 swap_em_state(&s_min, &s_b);
1525 /* Print it if necessary */
1528 if (mdrunOptions.verbose)
1530 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1531 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step,
1532 s_min->epot, s_min->fnorm / sqrtNumAtoms, s_min->fmax, s_min->a_fmax + 1);
1535 /* Store the new (lower) energies */
1536 matrix nullBox = {};
1537 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1538 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1539 nullptr, vir, pres, nullptr, mu_tot, constr);
1541 do_log = do_per_step(step, inputrec->nstlog);
1542 do_ene = do_per_step(step, inputrec->nstenergy);
1544 imdSession->fillEnergyRecord(step, TRUE);
1548 EnergyOutput::printHeader(fplog, step, step);
1550 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1551 do_log ? fplog : nullptr, step, step,
1552 fr->fcdata.get(), nullptr);
1555 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1556 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0))
1558 imdSession->sendPositionsAndEnergies();
1561 /* Stop when the maximum force lies below tolerance.
1562 * If we have reached machine precision, converged is already set to true.
1564 converged = converged || (s_min->fmax < inputrec->em_tol);
1566 } /* End of the loop */
1570 step--; /* we never took that last step in this case */
1572 if (s_min->fmax > inputrec->em_tol)
1576 warn_step(fplog, inputrec->em_tol, s_min->fmax, step - 1 == number_steps, FALSE);
1583 /* If we printed energy and/or logfile last step (which was the last step)
1584 * we don't have to do it again, but otherwise print the final values.
1588 /* Write final value to log since we didn't do anything the last step */
1589 EnergyOutput::printHeader(fplog, step, step);
1591 if (!do_ene || !do_log)
1593 /* Write final energy file entries */
1594 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1595 !do_log ? fplog : nullptr, step, step,
1596 fr->fcdata.get(), nullptr);
1600 /* Print some stuff... */
1603 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1607 * For accurate normal mode calculation it is imperative that we
1608 * store the last conformation into the full precision binary trajectory.
1610 * However, we should only do it if we did NOT already write this step
1611 * above (which we did if do_x or do_f was true).
1613 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1614 * in the trajectory with the same step number.
1616 do_x = !do_per_step(step, inputrec->nstxout);
1617 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1619 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec,
1620 step, s_min, state_global, observablesHistory);
1625 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1626 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1627 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1629 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1632 finish_em(cr, outf, walltime_accounting, wcycle);
1634 /* To print the actual number of steps we needed somewhere */
1635 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1639 void LegacySimulator::do_lbfgs()
1641 static const char* LBFGS = "Low-Memory BFGS Minimizer";
1643 gmx_localtop_t top(top_global->ffparams);
1644 gmx_global_stat_t gstat;
1645 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1646 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1647 real * rho, *alpha, *p, *s, **dx, **dg;
1648 real a, b, c, maxdelta, delta;
1650 real dgdx, dgdg, sq, yr, beta;
1652 rvec mu_tot = { 0 };
1653 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1655 int start, end, number_steps;
1656 int i, k, m, n, gf, step;
1658 auto mdatoms = mdAtoms->mdatoms();
1663 "Note that activating L-BFGS energy minimization via the "
1664 "integrator .mdp option and the command gmx mdrun may "
1665 "be available in a different form in a future version of GROMACS, "
1666 "e.g. gmx minimize and an .mdp option.");
1670 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1673 if (nullptr != constr)
1677 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1678 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1681 n = 3 * state_global->natoms;
1682 nmaxcorr = inputrec->nbfgscorr;
1687 snew(rho, nmaxcorr);
1688 snew(alpha, nmaxcorr);
1691 for (i = 0; i < nmaxcorr; i++)
1697 for (i = 0; i < nmaxcorr; i++)
1706 init_em(fplog, mdlog, LBFGS, cr, inputrec, imdSession, pull_work, state_global, top_global,
1707 &ems, &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, nullptr);
1708 const bool simulationsShareState = false;
1709 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
1710 mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1711 StartingBehavior::NewSimulation, simulationsShareState, ms);
1712 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
1713 false, StartingBehavior::NewSimulation, mdModulesNotifier);
1716 end = mdatoms->homenr;
1718 /* We need 4 working states */
1719 em_state_t s0{}, s1{}, s2{}, s3{};
1720 em_state_t* sa = &s0;
1721 em_state_t* sb = &s1;
1722 em_state_t* sc = &s2;
1723 em_state_t* last = &s3;
1724 /* Initialize by copying the state from ems (we could skip x and f here) */
1729 /* Print to log file */
1730 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1732 do_log = do_ene = do_x = do_f = TRUE;
1734 /* Max number of steps */
1735 number_steps = inputrec->nsteps;
1737 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1739 for (i = start; i < end; i++)
1741 if (mdatoms->cFREEZE)
1743 gf = mdatoms->cFREEZE[i];
1745 for (m = 0; m < DIM; m++)
1747 frozen[3 * i + m] = (inputrec->opts.nFreeze[gf][m] != 0);
1752 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1756 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1761 vsite->construct(state_global->x, 1, {}, state_global->box);
1764 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1765 /* do_force always puts the charge groups in the box and shifts again
1766 * We do not unshift, so molecules are always whole
1769 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
1770 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1771 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
1772 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1776 /* Copy stuff to the energy bin for easy printing etc. */
1777 matrix nullBox = {};
1778 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1779 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1780 nullptr, vir, pres, nullptr, mu_tot, constr);
1782 EnergyOutput::printHeader(fplog, step, step);
1783 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step,
1784 step, fr->fcdata.get(), nullptr);
1787 /* Set the initial step.
1788 * since it will be multiplied by the non-normalized search direction
1789 * vector (force vector the first time), we scale it by the
1790 * norm of the force.
1795 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1796 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1797 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1798 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
1799 fprintf(stderr, "\n");
1800 /* and copy to the log file too... */
1801 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1802 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1803 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
1804 fprintf(fplog, "\n");
1807 // Point is an index to the memory of search directions, where 0 is the first one.
1810 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1811 real* fInit = static_cast<real*>(ems.f.view().force().data()[0]);
1812 for (i = 0; i < n; i++)
1816 dx[point][i] = fInit[i]; /* Initial search direction */
1824 // Stepsize will be modified during the search, and actually it is not critical
1825 // (the main efficiency in the algorithm comes from changing directions), but
1826 // we still need an initial value, so estimate it as the inverse of the norm
1827 // so we take small steps where the potential fluctuates a lot.
1828 stepsize = 1.0 / ems.fnorm;
1830 /* Start the loop over BFGS steps.
1831 * Each successful step is counted, and we continue until
1832 * we either converge or reach the max number of steps.
1837 /* Set the gradient from the force */
1839 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1842 /* Write coordinates if necessary */
1843 do_x = do_per_step(step, inputrec->nstxout);
1844 do_f = do_per_step(step, inputrec->nstfout);
1849 mdof_flags |= MDOF_X;
1854 mdof_flags |= MDOF_F;
1859 mdof_flags |= MDOF_IMD;
1862 gmx::WriteCheckpointDataHolder checkpointDataHolder;
1863 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global->natoms, step,
1864 static_cast<real>(step), &ems.s, state_global, observablesHistory,
1865 ems.f.view().force(), &checkpointDataHolder);
1867 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1869 /* make s a pointer to current search direction - point=0 first time we get here */
1872 real* xx = static_cast<real*>(ems.s.x.rvec_array()[0]);
1873 real* ff = static_cast<real*>(ems.f.view().force().data()[0]);
1875 // calculate line gradient in position A
1876 for (gpa = 0, i = 0; i < n; i++)
1878 gpa -= s[i] * ff[i];
1881 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1882 * relative change in coordinate is smaller than precision
1884 for (minstep = 0, i = 0; i < n; i++)
1892 minstep += tmp * tmp;
1894 minstep = GMX_REAL_EPS / sqrt(minstep / n);
1896 if (stepsize < minstep)
1902 // Before taking any steps along the line, store the old position
1904 real* lastx = static_cast<real*>(last->s.x.data()[0]);
1905 real* lastf = static_cast<real*>(last->f.view().force().data()[0]);
1910 /* Take a step downhill.
1911 * In theory, we should find the actual minimum of the function in this
1912 * direction, somewhere along the line.
1913 * That is quite possible, but it turns out to take 5-10 function evaluations
1914 * for each line. However, we dont really need to find the exact minimum -
1915 * it is much better to start a new BFGS step in a modified direction as soon
1916 * as we are close to it. This will save a lot of energy evaluations.
1918 * In practice, we just try to take a single step.
1919 * If it worked (i.e. lowered the energy), we increase the stepsize but
1920 * continue straight to the next BFGS step without trying to find any minimum,
1921 * i.e. we change the search direction too. If the line was smooth, it is
1922 * likely we are in a smooth region, and then it makes sense to take longer
1923 * steps in the modified search direction too.
1925 * If it didn't work (higher energy), there must be a minimum somewhere between
1926 * the old position and the new one. Then we need to start by finding a lower
1927 * value before we change search direction. Since the energy was apparently
1928 * quite rough, we need to decrease the step size.
1930 * Due to the finite numerical accuracy, it turns out that it is a good idea
1931 * to accept a SMALL increase in energy, if the derivative is still downhill.
1932 * This leads to lower final energies in the tests I've done. / Erik
1935 // State "A" is the first position along the line.
1936 // reference position along line is initially zero
1939 // Check stepsize first. We do not allow displacements
1940 // larger than emstep.
1944 // Pick a new position C by adding stepsize to A.
1947 // Calculate what the largest change in any individual coordinate
1948 // would be (translation along line * gradient along line)
1950 for (i = 0; i < n; i++)
1953 if (delta > maxdelta)
1958 // If any displacement is larger than the stepsize limit, reduce the step
1959 if (maxdelta > inputrec->em_stepsize)
1963 } while (maxdelta > inputrec->em_stepsize);
1965 // Take a trial step and move the coordinate array xc[] to position C
1966 real* xc = static_cast<real*>(sc->s.x.rvec_array()[0]);
1967 for (i = 0; i < n; i++)
1969 xc[i] = lastx[i] + c * s[i];
1973 // Calculate energy for the trial step in position C
1974 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
1976 // Calc line gradient in position C
1977 real* fc = static_cast<real*>(sc->f.view().force()[0]);
1978 for (gpc = 0, i = 0; i < n; i++)
1980 gpc -= s[i] * fc[i]; /* f is negative gradient, thus the sign */
1982 /* Sum the gradient along the line across CPUs */
1985 gmx_sumd(1, &gpc, cr);
1988 // This is the max amount of increase in energy we tolerate.
1989 // By allowing VERY small changes (close to numerical precision) we
1990 // frequently find even better (lower) final energies.
1991 tmp = std::sqrt(GMX_REAL_EPS) * fabs(sa->epot);
1993 // Accept the step if the energy is lower in the new position C (compared to A),
1994 // or if it is not significantly higher and the line derivative is still negative.
1995 foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
1996 // If true, great, we found a better energy. We no longer try to alter the
1997 // stepsize, but simply accept this new better position. The we select a new
1998 // search direction instead, which will be much more efficient than continuing
1999 // to take smaller steps along a line. Set fnorm based on the new C position,
2000 // which will be used to update the stepsize to 1/fnorm further down.
2002 // If false, the energy is NOT lower in point C, i.e. it will be the same
2003 // or higher than in point A. In this case it is pointless to move to point C,
2004 // so we will have to do more iterations along the same line to find a smaller
2005 // value in the interval [A=0.0,C].
2006 // Here, A is still 0.0, but that will change when we do a search in the interval
2007 // [0.0,C] below. That search we will do by interpolation or bisection rather
2008 // than with the stepsize, so no need to modify it. For the next search direction
2009 // it will be reset to 1/fnorm anyway.
2013 // OK, if we didn't find a lower value we will have to locate one now - there must
2014 // be one in the interval [a,c].
2015 // The same thing is valid here, though: Don't spend dozens of iterations to find
2016 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2017 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2018 // I also have a safeguard for potentially really pathological functions so we never
2019 // take more than 20 steps before we give up.
2020 // If we already found a lower value we just skip this step and continue to the update.
2025 // Select a new trial point B in the interval [A,C].
2026 // If the derivatives at points a & c have different sign we interpolate to zero,
2027 // otherwise just do a bisection since there might be multiple minima/maxima
2028 // inside the interval.
2029 if (gpa < 0 && gpc > 0)
2031 b = a + gpa * (a - c) / (gpc - gpa);
2038 /* safeguard if interpolation close to machine accuracy causes errors:
2039 * never go outside the interval
2041 if (b <= a || b >= c)
2046 // Take a trial step to point B
2047 real* xb = static_cast<real*>(sb->s.x.rvec_array()[0]);
2048 for (i = 0; i < n; i++)
2050 xb[i] = lastx[i] + b * s[i];
2054 // Calculate energy for the trial step in point B
2055 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2058 // Calculate gradient in point B
2059 real* fb = static_cast<real*>(sb->f.view().force()[0]);
2060 for (gpb = 0, i = 0; i < n; i++)
2062 gpb -= s[i] * fb[i]; /* f is negative gradient, thus the sign */
2064 /* Sum the gradient along the line across CPUs */
2067 gmx_sumd(1, &gpb, cr);
2070 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2071 // at the new point B, and rename the endpoints of this new interval A and C.
2074 /* Replace c endpoint with b */
2076 /* copy state b to c */
2081 /* Replace a endpoint with b */
2083 /* copy state b to a */
2088 * Stop search as soon as we find a value smaller than the endpoints,
2089 * or if the tolerance is below machine precision.
2090 * Never run more than 20 steps, no matter what.
2093 } while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2095 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2097 /* OK. We couldn't find a significantly lower energy.
2098 * If ncorr==0 this was steepest descent, and then we give up.
2099 * If not, reset memory to restart as steepest descent before quitting.
2111 /* Search in gradient direction */
2112 for (i = 0; i < n; i++)
2114 dx[point][i] = ff[i];
2116 /* Reset stepsize */
2117 stepsize = 1.0 / fnorm;
2122 /* Select min energy state of A & C, put the best in xx/ff/Epot
2124 if (sc->epot < sa->epot)
2145 /* Update the memory information, and calculate a new
2146 * approximation of the inverse hessian
2149 /* Have new data in Epot, xx, ff */
2150 if (ncorr < nmaxcorr)
2155 for (i = 0; i < n; i++)
2157 dg[point][i] = lastf[i] - ff[i];
2158 dx[point][i] *= step_taken;
2163 for (i = 0; i < n; i++)
2165 dgdg += dg[point][i] * dg[point][i];
2166 dgdx += dg[point][i] * dx[point][i];
2171 rho[point] = 1.0 / dgdx;
2174 if (point >= nmaxcorr)
2180 for (i = 0; i < n; i++)
2187 /* Recursive update. First go back over the memory points */
2188 for (k = 0; k < ncorr; k++)
2197 for (i = 0; i < n; i++)
2199 sq += dx[cp][i] * p[i];
2202 alpha[cp] = rho[cp] * sq;
2204 for (i = 0; i < n; i++)
2206 p[i] -= alpha[cp] * dg[cp][i];
2210 for (i = 0; i < n; i++)
2215 /* And then go forward again */
2216 for (k = 0; k < ncorr; k++)
2219 for (i = 0; i < n; i++)
2221 yr += p[i] * dg[cp][i];
2224 beta = rho[cp] * yr;
2225 beta = alpha[cp] - beta;
2227 for (i = 0; i < n; i++)
2229 p[i] += beta * dx[cp][i];
2239 for (i = 0; i < n; i++)
2243 dx[point][i] = p[i];
2251 /* Print it if necessary */
2254 if (mdrunOptions.verbose)
2256 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2257 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step,
2258 ems.epot, ems.fnorm / sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2261 /* Store the new (lower) energies */
2262 matrix nullBox = {};
2263 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
2264 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
2265 nullptr, vir, pres, nullptr, mu_tot, constr);
2267 do_log = do_per_step(step, inputrec->nstlog);
2268 do_ene = do_per_step(step, inputrec->nstenergy);
2270 imdSession->fillEnergyRecord(step, TRUE);
2274 EnergyOutput::printHeader(fplog, step, step);
2276 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2277 do_log ? fplog : nullptr, step, step,
2278 fr->fcdata.get(), nullptr);
2281 /* Send x and E to IMD client, if bIMD is TRUE. */
2282 if (imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0) && MASTER(cr))
2284 imdSession->sendPositionsAndEnergies();
2287 // Reset stepsize in we are doing more iterations
2290 /* Stop when the maximum force lies below tolerance.
2291 * If we have reached machine precision, converged is already set to true.
2293 converged = converged || (ems.fmax < inputrec->em_tol);
2295 } /* End of the loop */
2299 step--; /* we never took that last step in this case */
2301 if (ems.fmax > inputrec->em_tol)
2305 warn_step(fplog, inputrec->em_tol, ems.fmax, step - 1 == number_steps, FALSE);
2310 /* If we printed energy and/or logfile last step (which was the last step)
2311 * we don't have to do it again, but otherwise print the final values.
2313 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2315 EnergyOutput::printHeader(fplog, step, step);
2317 if (!do_ene || !do_log) /* Write final energy file entries */
2319 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2320 !do_log ? fplog : nullptr, step, step, fr->fcdata.get(),
2324 /* Print some stuff... */
2327 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2331 * For accurate normal mode calculation it is imperative that we
2332 * store the last conformation into the full precision binary trajectory.
2334 * However, we should only do it if we did NOT already write this step
2335 * above (which we did if do_x or do_f was true).
2337 do_x = !do_per_step(step, inputrec->nstxout);
2338 do_f = !do_per_step(step, inputrec->nstfout);
2339 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec,
2340 step, &ems, state_global, observablesHistory);
2344 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2345 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2346 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2348 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2351 finish_em(cr, outf, walltime_accounting, wcycle);
2353 /* To print the actual number of steps we needed somewhere */
2354 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2357 void LegacySimulator::do_steep()
2359 const char* SD = "Steepest Descents";
2360 gmx_localtop_t top(top_global->ffparams);
2361 gmx_global_stat_t gstat;
2364 gmx_bool bDone, bAbort, do_x, do_f;
2366 rvec mu_tot = { 0 };
2369 int steps_accepted = 0;
2370 auto mdatoms = mdAtoms->mdatoms();
2375 "Note that activating steepest-descent energy minimization via the "
2376 "integrator .mdp option and the command gmx mdrun may "
2377 "be available in a different form in a future version of GROMACS, "
2378 "e.g. gmx minimize and an .mdp option.");
2380 /* Create 2 states on the stack and extract pointers that we will swap */
2381 em_state_t s0{}, s1{};
2382 em_state_t* s_min = &s0;
2383 em_state_t* s_try = &s1;
2385 /* Init em and store the local state in s_try */
2386 init_em(fplog, mdlog, SD, cr, inputrec, imdSession, pull_work, state_global, top_global, s_try,
2387 &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, nullptr);
2388 const bool simulationsShareState = false;
2389 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
2390 mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2391 StartingBehavior::NewSimulation, simulationsShareState, ms);
2392 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
2393 false, StartingBehavior::NewSimulation, mdModulesNotifier);
2395 /* Print to log file */
2396 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2398 /* Set variables for stepsize (in nm). This is the largest
2399 * step that we are going to make in any direction.
2401 ustep = inputrec->em_stepsize;
2404 /* Max number of steps */
2405 nsteps = inputrec->nsteps;
2409 /* Print to the screen */
2410 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2414 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2416 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2417 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2418 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2420 /**** HERE STARTS THE LOOP ****
2421 * count is the counter for the number of steps
2422 * bDone will be TRUE when the minimization has converged
2423 * bAbort will be TRUE when nsteps steps have been performed or when
2424 * the stepsize becomes smaller than is reasonable for machine precision
2429 while (!bDone && !bAbort)
2431 bAbort = (nsteps >= 0) && (count == nsteps);
2433 /* set new coordinates, except for first step */
2434 bool validStep = true;
2437 validStep = do_em_step(cr, inputrec, mdatoms, s_min, stepsize,
2438 s_min->f.view().forceWithPadding(), s_try, constr, count);
2443 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2447 // Signal constraint error during stepping with energy=inf
2448 s_try->epot = std::numeric_limits<real>::infinity();
2453 EnergyOutput::printHeader(fplog, count, count);
2458 s_min->epot = s_try->epot;
2461 /* Print it if necessary */
2464 if (mdrunOptions.verbose)
2466 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2467 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax + 1,
2468 ((count == 0) || (s_try->epot < s_min->epot)) ? '\n' : '\r');
2472 if ((count == 0) || (s_try->epot < s_min->epot))
2474 /* Store the new (lower) energies */
2475 matrix nullBox = {};
2476 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(count), mdatoms->tmass,
2477 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
2478 nullptr, vir, pres, nullptr, mu_tot, constr);
2480 imdSession->fillEnergyRecord(count, TRUE);
2482 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2483 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2484 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, do_dr, do_or,
2485 fplog, count, count, fr->fcdata.get(), nullptr);
2490 /* Now if the new energy is smaller than the previous...
2491 * or if this is the first step!
2492 * or if we did random steps!
2495 if ((count == 0) || (s_try->epot < s_min->epot))
2499 /* Test whether the convergence criterion is met... */
2500 bDone = (s_try->fmax < inputrec->em_tol);
2502 /* Copy the arrays for force, positions and energy */
2503 /* The 'Min' array always holds the coords and forces of the minimal
2505 swap_em_state(&s_min, &s_try);
2511 /* Write to trn, if necessary */
2512 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2513 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2514 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, count, s_min,
2515 state_global, observablesHistory);
2519 /* If energy is not smaller make the step smaller... */
2522 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2524 /* Reload the old state */
2525 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession,
2526 pull_work, s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
2530 // If the force is very small after finishing minimization,
2531 // we risk dividing by zero when calculating the step size.
2532 // So we check first if the minimization has stopped before
2533 // trying to obtain a new step size.
2536 /* Determine new step */
2537 stepsize = ustep / s_min->fmax;
2540 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2542 if (count == nsteps || ustep < 1e-12)
2544 if (count == nsteps || ustep < 1e-6)
2549 warn_step(fplog, inputrec->em_tol, s_min->fmax, count == nsteps, constr != nullptr);
2554 /* Send IMD energies and positions, if bIMD is TRUE. */
2555 if (imdSession->run(count, TRUE, MASTER(cr) ? state_global->box : nullptr,
2556 MASTER(cr) ? state_global->x.rvec_array() : nullptr, 0)
2559 imdSession->sendPositionsAndEnergies();
2563 } /* End of the loop */
2565 /* Print some data... */
2568 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2570 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2571 top_global, inputrec, count, s_min, state_global, observablesHistory);
2575 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2577 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2578 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2581 finish_em(cr, outf, walltime_accounting, wcycle);
2583 /* To print the actual number of steps we needed somewhere */
2584 inputrec->nsteps = count;
2586 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2589 void LegacySimulator::do_nm()
2591 const char* NM = "Normal Mode Analysis";
2593 gmx_localtop_t top(top_global->ffparams);
2594 gmx_global_stat_t gstat;
2596 rvec mu_tot = { 0 };
2598 gmx_bool bSparse; /* use sparse matrix storage format */
2600 gmx_sparsematrix_t* sparse_matrix = nullptr;
2601 real* full_matrix = nullptr;
2603 /* added with respect to mdrun */
2605 real der_range = 10.0 * std::sqrt(GMX_REAL_EPS);
2607 bool bIsMaster = MASTER(cr);
2608 auto mdatoms = mdAtoms->mdatoms();
2613 "Note that activating normal-mode analysis via the integrator "
2614 ".mdp option and the command gmx mdrun may "
2615 "be available in a different form in a future version of GROMACS, "
2616 "e.g. gmx normal-modes.");
2618 if (constr != nullptr)
2622 "Constraints present with Normal Mode Analysis, this combination is not supported");
2625 gmx_shellfc_t* shellfc;
2627 em_state_t state_work{};
2629 /* Init em and store the local state in state_minimum */
2630 init_em(fplog, mdlog, NM, cr, inputrec, imdSession, pull_work, state_global, top_global,
2631 &state_work, &top, nrnb, fr, mdAtoms, &gstat, vsite, constr, &shellfc);
2632 const bool simulationsShareState = false;
2633 gmx_mdoutf* outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider,
2634 mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2635 StartingBehavior::NewSimulation, simulationsShareState, ms);
2637 std::vector<int> atom_index = get_atom_index(top_global);
2638 std::vector<gmx::RVec> fneg(atom_index.size(), { 0, 0, 0 });
2639 snew(dfdx, atom_index.size());
2645 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2646 " which MIGHT not be accurate enough for normal mode analysis.\n"
2647 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2648 " are fairly modest even if you recompile in double precision.\n\n");
2652 /* Check if we can/should use sparse storage format.
2654 * Sparse format is only useful when the Hessian itself is sparse, which it
2655 * will be when we use a cutoff.
2656 * For small systems (n<1000) it is easier to always use full matrix format, though.
2658 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2660 GMX_LOG(mdlog.warning)
2661 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2664 else if (atom_index.size() < 1000)
2666 GMX_LOG(mdlog.warning)
2667 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2673 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2677 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2678 sz = DIM * atom_index.size();
2680 fprintf(stderr, "Allocating Hessian memory...\n\n");
2684 sparse_matrix = gmx_sparsematrix_init(sz);
2685 sparse_matrix->compressed_symmetric = TRUE;
2689 snew(full_matrix, sz * sz);
2692 /* Write start time and temperature */
2693 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2695 /* fudge nr of steps to nr of atoms */
2696 inputrec->nsteps = atom_index.size() * 2;
2700 fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2701 *(top_global->name), inputrec->nsteps);
2704 nnodes = cr->nnodes;
2706 /* Make evaluate_energy do a single node force calculation */
2708 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2709 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2710 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2711 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2712 cr->nnodes = nnodes;
2714 /* if forces are not small, warn user */
2715 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2717 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2718 if (state_work.fmax > 1.0e-3)
2720 GMX_LOG(mdlog.warning)
2722 "The force is probably not small enough to "
2723 "ensure that you are at a minimum.\n"
2724 "Be aware that negative eigenvalues may occur\n"
2725 "when the resulting matrix is diagonalized.");
2728 /***********************************************************
2730 * Loop over all pairs in matrix
2732 * do_force called twice. Once with positive and
2733 * once with negative displacement
2735 ************************************************************/
2737 /* Steps are divided one by one over the nodes */
2739 auto state_work_x = makeArrayRef(state_work.s.x);
2740 auto state_work_f = state_work.f.view().force();
2741 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
2743 size_t atom = atom_index[aid];
2744 for (size_t d = 0; d < DIM; d++)
2747 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2750 x_min = state_work_x[atom][d];
2752 for (unsigned int dx = 0; (dx < 2); dx++)
2756 state_work_x[atom][d] = x_min - der_range;
2760 state_work_x[atom][d] = x_min + der_range;
2763 /* Make evaluate_energy do a single node force calculation */
2767 /* Now is the time to relax the shells */
2768 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, nullptr, step, inputrec,
2769 imdSession, pull_work, bNS, force_flags, &top, constr, enerd,
2770 state_work.s.natoms, state_work.s.x.arrayRefWithPadding(),
2771 state_work.s.v.arrayRefWithPadding(), state_work.s.box,
2772 state_work.s.lambda, &state_work.s.hist, &state_work.f.view(),
2773 vir, mdatoms, nrnb, wcycle, shellfc, fr, runScheduleWork, t,
2774 mu_tot, vsite, DDBalanceRegionHandler(nullptr));
2780 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid * 2 + dx, FALSE);
2783 cr->nnodes = nnodes;
2787 std::copy(state_work_f.begin(), state_work_f.begin() + atom_index.size(),
2792 /* x is restored to original */
2793 state_work_x[atom][d] = x_min;
2795 for (size_t j = 0; j < atom_index.size(); j++)
2797 for (size_t k = 0; (k < DIM); k++)
2799 dfdx[j][k] = -(state_work_f[atom_index[j]][k] - fneg[j][k]) / (2 * der_range);
2806 # define mpi_type GMX_MPI_REAL
2807 MPI_Send(dfdx[0], atom_index.size() * DIM, mpi_type, MASTER(cr), cr->nodeid,
2808 cr->mpi_comm_mygroup);
2813 for (index node = 0; (node < nnodes && aid + node < ssize(atom_index)); node++)
2819 MPI_Recv(dfdx[0], atom_index.size() * DIM, mpi_type, node, node,
2820 cr->mpi_comm_mygroup, &stat);
2825 row = (aid + node) * DIM + d;
2827 for (size_t j = 0; j < atom_index.size(); j++)
2829 for (size_t k = 0; k < DIM; k++)
2835 if (col >= row && dfdx[j][k] != 0.0)
2837 gmx_sparsematrix_increment_value(sparse_matrix, row, col, dfdx[j][k]);
2842 full_matrix[row * sz + col] = dfdx[j][k];
2849 if (mdrunOptions.verbose && fplog)
2854 /* write progress */
2855 if (bIsMaster && mdrunOptions.verbose)
2857 fprintf(stderr, "\rFinished step %d out of %td",
2858 std::min<int>(atom + nnodes, atom_index.size()), ssize(atom_index));
2865 fprintf(stderr, "\n\nWriting Hessian...\n");
2866 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2869 finish_em(cr, outf, walltime_accounting, wcycle);
2871 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size() * 2);