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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file defines integrators for energy minimization
41 * \author Berk Hess <hess@kth.se>
42 * \author Erik Lindahl <erik@kth.se>
43 * \ingroup module_mdrun
56 #include "gromacs/commandline/filenm.h"
57 #include "gromacs/domdec/collect.h"
58 #include "gromacs/domdec/dlbtiming.h"
59 #include "gromacs/domdec/domdec.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/domdec/mdsetup.h"
62 #include "gromacs/domdec/partition.h"
63 #include "gromacs/ewald/pme.h"
64 #include "gromacs/fileio/confio.h"
65 #include "gromacs/fileio/mtxio.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gmxlib/nrnb.h"
68 #include "gromacs/imd/imd.h"
69 #include "gromacs/linearalgebra/sparsematrix.h"
70 #include "gromacs/listed_forces/manage_threading.h"
71 #include "gromacs/math/functions.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/mdlib/constr.h"
74 #include "gromacs/mdlib/dispersioncorrection.h"
75 #include "gromacs/mdlib/ebin.h"
76 #include "gromacs/mdlib/enerdata_utils.h"
77 #include "gromacs/mdlib/energyoutput.h"
78 #include "gromacs/mdlib/force.h"
79 #include "gromacs/mdlib/forcerec.h"
80 #include "gromacs/mdlib/gmx_omp_nthreads.h"
81 #include "gromacs/mdlib/md_support.h"
82 #include "gromacs/mdlib/mdatoms.h"
83 #include "gromacs/mdlib/stat.h"
84 #include "gromacs/mdlib/tgroup.h"
85 #include "gromacs/mdlib/trajectory_writing.h"
86 #include "gromacs/mdlib/update.h"
87 #include "gromacs/mdlib/vsite.h"
88 #include "gromacs/mdrunutility/handlerestart.h"
89 #include "gromacs/mdrunutility/printtime.h"
90 #include "gromacs/mdtypes/commrec.h"
91 #include "gromacs/mdtypes/inputrec.h"
92 #include "gromacs/mdtypes/md_enums.h"
93 #include "gromacs/mdtypes/mdrunoptions.h"
94 #include "gromacs/mdtypes/state.h"
95 #include "gromacs/pbcutil/mshift.h"
96 #include "gromacs/pbcutil/pbc.h"
97 #include "gromacs/timing/wallcycle.h"
98 #include "gromacs/timing/walltime_accounting.h"
99 #include "gromacs/topology/mtop_util.h"
100 #include "gromacs/topology/topology.h"
101 #include "gromacs/utility/cstringutil.h"
102 #include "gromacs/utility/exceptions.h"
103 #include "gromacs/utility/fatalerror.h"
104 #include "gromacs/utility/logger.h"
105 #include "gromacs/utility/smalloc.h"
107 #include "legacysimulator.h"
110 using gmx::MdrunScheduleWorkload;
112 //! Utility structure for manipulating states during EM
115 //! Copy of the global state
118 PaddedHostVector<gmx::RVec> f;
121 //! Norm of the force
129 //! Print the EM starting conditions
130 static void print_em_start(FILE* fplog,
132 gmx_walltime_accounting_t walltime_accounting,
133 gmx_wallcycle_t wcycle,
136 walltime_accounting_start_time(walltime_accounting);
137 wallcycle_start(wcycle, ewcRUN);
138 print_start(fplog, cr, walltime_accounting, name);
141 //! Stop counting time for EM
142 static void em_time_end(gmx_walltime_accounting_t walltime_accounting, gmx_wallcycle_t wcycle)
144 wallcycle_stop(wcycle, ewcRUN);
146 walltime_accounting_end_time(walltime_accounting);
149 //! Printing a log file and console header
150 static void sp_header(FILE* out, const char* minimizer, real ftol, int nsteps)
153 fprintf(out, "%s:\n", minimizer);
154 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
155 fprintf(out, " Number of steps = %12d\n", nsteps);
158 //! Print warning message
159 static void warn_step(FILE* fp, real ftol, real fmax, gmx_bool bLastStep, gmx_bool bConstrain)
161 constexpr bool realIsDouble = GMX_DOUBLE;
164 if (!std::isfinite(fmax))
167 "\nEnergy minimization has stopped because the force "
168 "on at least one atom is not finite. This usually means "
169 "atoms are overlapping. Modify the input coordinates to "
170 "remove atom overlap or use soft-core potentials with "
171 "the free energy code to avoid infinite forces.\n%s",
172 !realIsDouble ? "You could also be lucky that switching to double precision "
173 "is sufficient to obtain finite forces.\n"
179 "\nEnergy minimization reached the maximum number "
180 "of steps before the forces reached the requested "
181 "precision Fmax < %g.\n",
187 "\nEnergy minimization has stopped, but the forces have "
188 "not converged to the requested precision Fmax < %g (which "
189 "may not be possible for your system). It stopped "
190 "because the algorithm tried to make a new step whose size "
191 "was too small, or there was no change in the energy since "
192 "last step. Either way, we regard the minimization as "
193 "converged to within the available machine precision, "
194 "given your starting configuration and EM parameters.\n%s%s",
196 !realIsDouble ? "\nDouble precision normally gives you higher accuracy, but "
197 "this is often not needed for preparing to run molecular "
200 bConstrain ? "You might need to increase your constraint accuracy, or turn\n"
201 "off constraints altogether (set constraints = none in mdp file)\n"
205 fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
206 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
209 //! Print message about convergence of the EM
210 static void print_converged(FILE* fp,
216 const em_state_t* ems,
219 char buf[STEPSTRSIZE];
223 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n", alg, ftol, gmx_step_str(count, buf));
225 else if (count < nsteps)
228 "\n%s converged to machine precision in %s steps,\n"
229 "but did not reach the requested Fmax < %g.\n",
230 alg, gmx_step_str(count, buf), ftol);
234 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n", alg, ftol,
235 gmx_step_str(count, buf));
239 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
240 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
241 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm / sqrtNumAtoms);
243 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
244 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
245 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm / sqrtNumAtoms);
249 //! Compute the norm and max of the force array in parallel
250 static void get_f_norm_max(const t_commrec* cr,
260 int la_max, a_max, start, end, i, m, gf;
262 /* This routine finds the largest force and returns it.
263 * On parallel machines the global max is taken.
269 end = mdatoms->homenr;
270 if (mdatoms->cFREEZE)
272 for (i = start; i < end; i++)
274 gf = mdatoms->cFREEZE[i];
276 for (m = 0; m < DIM; m++)
278 if (!opts->nFreeze[gf][m])
280 fam += gmx::square(f[i][m]);
293 for (i = start; i < end; i++)
305 if (la_max >= 0 && DOMAINDECOMP(cr))
307 a_max = cr->dd->globalAtomIndices[la_max];
315 snew(sum, 2 * cr->nnodes + 1);
316 sum[2 * cr->nodeid] = fmax2;
317 sum[2 * cr->nodeid + 1] = a_max;
318 sum[2 * cr->nnodes] = fnorm2;
319 gmx_sumd(2 * cr->nnodes + 1, sum, cr);
320 fnorm2 = sum[2 * cr->nnodes];
321 /* Determine the global maximum */
322 for (i = 0; i < cr->nnodes; i++)
324 if (sum[2 * i] > fmax2)
327 a_max = gmx::roundToInt(sum[2 * i + 1]);
335 *fnorm = sqrt(fnorm2);
347 //! Compute the norm of the force
348 static void get_state_f_norm_max(const t_commrec* cr, t_grpopts* opts, t_mdatoms* mdatoms, em_state_t* ems)
350 get_f_norm_max(cr, opts, mdatoms, ems->f.rvec_array(), &ems->fnorm, &ems->fmax, &ems->a_fmax);
353 //! Initialize the energy minimization
354 static void init_em(FILE* fplog,
355 const gmx::MDLogger& mdlog,
359 gmx::ImdSession* imdSession,
361 t_state* state_global,
362 gmx_mtop_t* top_global,
368 gmx::MDAtoms* mdAtoms,
369 gmx_global_stat_t* gstat,
371 gmx::Constraints* constr,
372 gmx_shellfc_t** shellfc)
378 fprintf(fplog, "Initiating %s\n", title);
383 state_global->ngtc = 0;
385 initialize_lambdas(fplog, *ir, MASTER(cr), &(state_global->fep_state), state_global->lambda, nullptr);
389 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
391 *shellfc = init_shell_flexcon(stdout, top_global, constr ? constr->numFlexibleConstraints() : 0,
392 ir->nstcalcenergy, DOMAINDECOMP(cr));
396 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI),
397 "This else currently only handles energy minimizers, consider if your algorithm "
398 "needs shell/flexible-constraint support");
400 /* With energy minimization, shells and flexible constraints are
401 * automatically minimized when treated like normal DOFS.
403 if (shellfc != nullptr)
409 auto mdatoms = mdAtoms->mdatoms();
410 if (DOMAINDECOMP(cr))
412 top->useInDomainDecomp_ = true;
413 dd_init_local_top(*top_global, top);
415 dd_init_local_state(cr->dd, state_global, &ems->s);
417 /* Distribute the charge groups over the nodes from the master node */
418 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1, state_global, *top_global, ir,
419 imdSession, pull_work, &ems->s, &ems->f, mdAtoms, top, fr, vsite,
420 constr, nrnb, nullptr, FALSE);
421 dd_store_state(cr->dd, &ems->s);
427 state_change_natoms(state_global, state_global->natoms);
428 /* Just copy the state */
429 ems->s = *state_global;
430 state_change_natoms(&ems->s, ems->s.natoms);
431 ems->f.resizeWithPadding(ems->s.natoms);
433 mdAlgorithmsSetupAtomData(cr, ir, *top_global, top, fr, graph, mdAtoms, constr, vsite,
434 shellfc ? *shellfc : nullptr);
438 set_vsite_top(vsite, top, mdatoms);
442 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
446 // TODO how should this cross-module support dependency be managed?
447 if (ir->eConstrAlg == econtSHAKE && gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
449 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
450 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
453 if (!ir->bContinuation)
455 /* Constrain the starting coordinates */
457 constr->apply(TRUE, TRUE, -1, 0, 1.0, ems->s.x.rvec_array(), ems->s.x.rvec_array(),
458 nullptr, ems->s.box, ems->s.lambda[efptFEP], &dvdl_constr, nullptr,
459 nullptr, gmx::ConstraintVariable::Positions);
465 *gstat = global_stat_init(ir);
472 calc_shifts(ems->s.box, fr->shift_vec);
475 //! Finalize the minimization
476 static void finish_em(const t_commrec* cr,
478 gmx_walltime_accounting_t walltime_accounting,
479 gmx_wallcycle_t wcycle)
481 if (!thisRankHasDuty(cr, DUTY_PME))
483 /* Tell the PME only node to finish */
484 gmx_pme_send_finish(cr);
489 em_time_end(walltime_accounting, wcycle);
492 //! Swap two different EM states during minimization
493 static void swap_em_state(em_state_t** ems1, em_state_t** ems2)
502 //! Save the EM trajectory
503 static void write_em_traj(FILE* fplog,
509 gmx_mtop_t* top_global,
513 t_state* state_global,
514 ObservablesHistory* observablesHistory)
520 mdof_flags |= MDOF_X;
524 mdof_flags |= MDOF_F;
527 /* If we want IMD output, set appropriate MDOF flag */
530 mdof_flags |= MDOF_IMD;
533 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global->natoms, step,
534 static_cast<double>(step), &state->s, state_global,
535 observablesHistory, state->f);
537 if (confout != nullptr)
539 if (DOMAINDECOMP(cr))
541 /* If bX=true, x was collected to state_global in the call above */
544 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
545 dd_collect_vec(cr->dd, &state->s, state->s.x, globalXRef);
550 /* Copy the local state pointer */
551 state_global = &state->s;
556 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
558 /* Make molecules whole only for confout writing */
559 do_pbc_mtop(ir->ePBC, state->s.box, top_global, state_global->x.rvec_array());
562 write_sto_conf_mtop(confout, *top_global->name, top_global,
563 state_global->x.rvec_array(), nullptr, ir->ePBC, state->s.box);
568 //! \brief Do one minimization step
570 // \returns true when the step succeeded, false when a constraint error occurred
571 static bool do_em_step(const t_commrec* cr,
576 const PaddedHostVector<gmx::RVec>* force,
578 gmx::Constraints* constr,
585 int nthreads gmx_unused;
587 bool validStep = true;
592 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
594 gmx_incons("state mismatch in do_em_step");
597 s2->flags = s1->flags;
599 if (s2->natoms != s1->natoms)
601 state_change_natoms(s2, s1->natoms);
602 ems2->f.resizeWithPadding(s2->natoms);
604 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
606 s2->cg_gl.resize(s1->cg_gl.size());
609 copy_mat(s1->box, s2->box);
610 /* Copy free energy state */
611 s2->lambda = s1->lambda;
612 copy_mat(s1->box, s2->box);
617 nthreads = gmx_omp_nthreads_get(emntUpdate);
618 #pragma omp parallel num_threads(nthreads)
620 const rvec* x1 = s1->x.rvec_array();
621 rvec* x2 = s2->x.rvec_array();
622 const rvec* f = force->rvec_array();
625 #pragma omp for schedule(static) nowait
626 for (int i = start; i < end; i++)
634 for (int m = 0; m < DIM; m++)
636 if (ir->opts.nFreeze[gf][m])
642 x2[i][m] = x1[i][m] + a * f[i][m];
646 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
649 if (s2->flags & (1 << estCGP))
651 /* Copy the CG p vector */
652 const rvec* p1 = s1->cg_p.rvec_array();
653 rvec* p2 = s2->cg_p.rvec_array();
654 #pragma omp for schedule(static) nowait
655 for (int i = start; i < end; i++)
657 // Trivial OpenMP block that does not throw
658 copy_rvec(p1[i], p2[i]);
662 if (DOMAINDECOMP(cr))
664 /* OpenMP does not supported unsigned loop variables */
665 #pragma omp for schedule(static) nowait
666 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
668 s2->cg_gl[i] = s1->cg_gl[i];
673 if (DOMAINDECOMP(cr))
675 s2->ddp_count = s1->ddp_count;
676 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
682 validStep = constr->apply(TRUE, TRUE, count, 0, 1.0, s1->x.rvec_array(), s2->x.rvec_array(),
683 nullptr, s2->box, s2->lambda[efptBONDED], &dvdl_constr, nullptr,
684 nullptr, gmx::ConstraintVariable::Positions);
688 /* This global reduction will affect performance at high
689 * parallelization, but we can not really avoid it.
690 * But usually EM is not run at high parallelization.
692 int reductionBuffer = static_cast<int>(!validStep);
693 gmx_sumi(1, &reductionBuffer, cr);
694 validStep = (reductionBuffer == 0);
697 // We should move this check to the different minimizers
698 if (!validStep && ir->eI != eiSteep)
701 "The coordinates could not be constrained. Minimizer '%s' can not handle "
702 "constraint failures, use minimizer '%s' before using '%s'.",
703 EI(ir->eI), EI(eiSteep), EI(ir->eI));
710 //! Prepare EM for using domain decomposition parallellization
711 static void em_dd_partition_system(FILE* fplog,
712 const gmx::MDLogger& mdlog,
715 gmx_mtop_t* top_global,
717 gmx::ImdSession* imdSession,
721 gmx::MDAtoms* mdAtoms,
724 gmx::Constraints* constr,
726 gmx_wallcycle_t wcycle)
728 /* Repartition the domain decomposition */
729 dd_partition_system(fplog, mdlog, step, cr, FALSE, 1, nullptr, *top_global, ir, imdSession, pull_work,
730 &ems->s, &ems->f, mdAtoms, top, fr, vsite, constr, nrnb, wcycle, FALSE);
731 dd_store_state(cr->dd, &ems->s);
737 /*! \brief Class to handle the work of setting and doing an energy evaluation.
739 * This class is a mere aggregate of parameters to pass to evaluate an
740 * energy, so that future changes to names and types of them consume
741 * less time when refactoring other code.
743 * Aggregate initialization is used, for which the chief risk is that
744 * if a member is added at the end and not all initializer lists are
745 * updated, then the member will be value initialized, which will
746 * typically mean initialization to zero.
748 * Use a braced initializer list to construct one of these. */
749 class EnergyEvaluator
752 /*! \brief Evaluates an energy on the state in \c ems.
754 * \todo In practice, the same objects mu_tot, vir, and pres
755 * are always passed to this function, so we would rather have
756 * them as data members. However, their C-array types are
757 * unsuited for aggregate initialization. When the types
758 * improve, the call signature of this method can be reduced.
760 void run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst);
761 //! Handles logging (deprecated).
764 const gmx::MDLogger& mdlog;
765 //! Handles communication.
767 //! Coordinates multi-simulations.
768 const gmx_multisim_t* ms;
769 //! Holds the simulation topology.
770 gmx_mtop_t* top_global;
771 //! Holds the domain topology.
773 //! User input options.
774 t_inputrec* inputrec;
775 //! The Interactive Molecular Dynamics session.
776 gmx::ImdSession* imdSession;
777 //! The pull work object.
779 //! Manages flop accounting.
781 //! Manages wall cycle accounting.
782 gmx_wallcycle_t wcycle;
783 //! Coordinates global reduction.
784 gmx_global_stat_t gstat;
785 //! Handles virtual sites.
787 //! Handles constraints.
788 gmx::Constraints* constr;
789 //! Handles strange things.
791 //! Molecular graph for SHAKE.
793 //! Per-atom data for this domain.
794 gmx::MDAtoms* mdAtoms;
795 //! Handles how to calculate the forces.
797 //! Schedule of force-calculation work each step for this task.
798 MdrunScheduleWorkload* runScheduleWork;
799 //! Stores the computed energies.
800 gmx_enerdata_t* enerd;
803 void EnergyEvaluator::run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst)
807 tensor force_vir, shake_vir, ekin;
811 /* Set the time to the initial time, the time does not change during EM */
812 t = inputrec->init_t;
814 if (bFirst || (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
816 /* This is the first state or an old state used before the last ns */
822 if (inputrec->nstlist > 0)
830 construct_vsites(vsite, ems->s.x.rvec_array(), 1, nullptr, top->idef.iparams, top->idef.il,
831 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
834 if (DOMAINDECOMP(cr) && bNS)
836 /* Repartition the domain decomposition */
837 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession, pull_work,
838 ems, top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
841 /* Calc force & energy on new trial position */
842 /* do_force always puts the charge groups in the box and shifts again
843 * We do not unshift, so molecules are always whole in congrad.c
845 do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession, pull_work, count, nrnb, wcycle,
846 top, ems->s.box, ems->s.x.arrayRefWithPadding(), &ems->s.hist,
847 ems->f.arrayRefWithPadding(), force_vir, mdAtoms->mdatoms(), enerd, fcd, ems->s.lambda,
848 graph, 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, mu_tot, inputrec, nullptr, nullptr, nullptr,
863 1, &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 */
890 rvec* f_rvec = ems->f.rvec_array();
891 constr->apply(FALSE, FALSE, count, 0, 1.0, ems->s.x.rvec_array(), f_rvec, f_rvec,
892 ems->s.box, ems->s.lambda[efptBONDED], &dvdl_constr, nullptr, &shake_vir,
893 gmx::ConstraintVariable::ForceDispl);
894 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
895 m_add(force_vir, shake_vir, vir);
899 copy_mat(force_vir, vir);
903 enerd->term[F_PRES] = calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
905 sum_dhdl(enerd, ems->s.lambda, *inputrec->fepvals);
907 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
909 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
915 //! Parallel utility summing energies and forces
916 static double reorder_partsum(const t_commrec* cr,
918 gmx_mtop_t* top_global,
924 fprintf(debug, "Doing reorder_partsum\n");
927 const rvec* fm = s_min->f.rvec_array();
928 const rvec* fb = s_b->f.rvec_array();
930 /* Collect fm in a global vector fmg.
931 * This conflicts with the spirit of domain decomposition,
932 * but to fully optimize this a much more complicated algorithm is required.
934 const int natoms = top_global->natoms;
938 gmx::ArrayRef<const int> indicesMin = s_min->s.cg_gl;
940 for (int a : indicesMin)
942 copy_rvec(fm[i], fmg[a]);
945 gmx_sum(top_global->natoms * 3, fmg[0], cr);
947 /* Now we will determine the part of the sum for the cgs in state s_b */
948 gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
953 gmx::ArrayRef<unsigned char> grpnrFREEZE =
954 top_global->groups.groupNumbers[SimulationAtomGroupType::Freeze];
955 for (int a : indicesB)
957 if (!grpnrFREEZE.empty())
961 for (int m = 0; m < DIM; m++)
963 if (!opts->nFreeze[gf][m])
965 partsum += (fb[i][m] - fmg[a][m]) * fb[i][m];
976 //! Print some stuff, like beta, whatever that means.
977 static real pr_beta(const t_commrec* cr,
980 gmx_mtop_t* top_global,
986 /* This is just the classical Polak-Ribiere calculation of beta;
987 * it looks a bit complicated since we take freeze groups into account,
988 * and might have to sum it in parallel runs.
991 if (!DOMAINDECOMP(cr)
992 || (s_min->s.ddp_count == cr->dd->ddp_count && s_b->s.ddp_count == cr->dd->ddp_count))
994 const rvec* fm = s_min->f.rvec_array();
995 const rvec* fb = s_b->f.rvec_array();
998 /* This part of code can be incorrect with DD,
999 * since the atom ordering in s_b and s_min might differ.
1001 for (int i = 0; i < mdatoms->homenr; i++)
1003 if (mdatoms->cFREEZE)
1005 gf = mdatoms->cFREEZE[i];
1007 for (int m = 0; m < DIM; m++)
1009 if (!opts->nFreeze[gf][m])
1011 sum += (fb[i][m] - fm[i][m]) * fb[i][m];
1018 /* We need to reorder cgs while summing */
1019 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1023 gmx_sumd(1, &sum, cr);
1026 return sum / gmx::square(s_min->fnorm);
1032 void LegacySimulator::do_cg()
1034 const char* CG = "Polak-Ribiere Conjugate Gradients";
1037 gmx_global_stat_t gstat;
1039 double tmp, minstep;
1041 real a, b, c, beta = 0.0;
1044 gmx_bool converged, foundlower;
1045 rvec mu_tot = { 0 };
1046 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1048 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1049 int m, step, nminstep;
1050 auto mdatoms = mdAtoms->mdatoms();
1055 "Note that activating conjugate gradient energy minimization via the "
1056 "integrator .mdp option and the command gmx mdrun may "
1057 "be available in a different form in a future version of GROMACS, "
1058 "e.g. gmx minimize and an .mdp option.");
1064 // In CG, the state is extended with a search direction
1065 state_global->flags |= (1 << estCGP);
1067 // Ensure the extra per-atom state array gets allocated
1068 state_change_natoms(state_global, state_global->natoms);
1070 // Initialize the search direction to zero
1071 for (RVec& cg_p : state_global->cg_p)
1077 /* Create 4 states on the stack and extract pointers that we will swap */
1078 em_state_t s0{}, s1{}, s2{}, s3{};
1079 em_state_t* s_min = &s0;
1080 em_state_t* s_a = &s1;
1081 em_state_t* s_b = &s2;
1082 em_state_t* s_c = &s3;
1084 /* Init em and store the local state in s_min */
1085 init_em(fplog, mdlog, CG, cr, inputrec, imdSession, pull_work, state_global, top_global, s_min,
1086 &top, nrnb, fr, &graph, mdAtoms, &gstat, vsite, constr, nullptr);
1088 init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier,
1089 inputrec, top_global, nullptr, wcycle, StartingBehavior::NewSimulation);
1090 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
1091 false, StartingBehavior::NewSimulation, mdModulesNotifier);
1093 /* Print to log file */
1094 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1096 /* Max number of steps */
1097 number_steps = inputrec->nsteps;
1101 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1105 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1108 EnergyEvaluator energyEvaluator{
1109 fplog, mdlog, cr, ms, top_global, &top, inputrec,
1110 imdSession, pull_work, nrnb, wcycle, gstat, vsite, constr,
1111 fcd, graph, mdAtoms, fr, runScheduleWork, enerd
1113 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1114 /* do_force always puts the charge groups in the box and shifts again
1115 * We do not unshift, so molecules are always whole in congrad.c
1117 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1121 /* Copy stuff to the energy bin for easy printing etc. */
1122 matrix nullBox = {};
1123 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1124 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1125 nullptr, vir, pres, nullptr, mu_tot, constr);
1127 energyOutput.printHeader(fplog, step, step);
1128 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step,
1129 step, fcd, nullptr);
1132 /* Estimate/guess the initial stepsize */
1133 stepsize = inputrec->em_stepsize / s_min->fnorm;
1137 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1138 fprintf(stderr, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1139 fprintf(stderr, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1140 fprintf(stderr, "\n");
1141 /* and copy to the log file too... */
1142 fprintf(fplog, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1143 fprintf(fplog, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1144 fprintf(fplog, "\n");
1146 /* Start the loop over CG steps.
1147 * Each successful step is counted, and we continue until
1148 * we either converge or reach the max number of steps.
1151 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1154 /* start taking steps in a new direction
1155 * First time we enter the routine, beta=0, and the direction is
1156 * simply the negative gradient.
1159 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1160 rvec* pm = s_min->s.cg_p.rvec_array();
1161 const rvec* sfm = s_min->f.rvec_array();
1164 for (int i = 0; i < mdatoms->homenr; i++)
1166 if (mdatoms->cFREEZE)
1168 gf = mdatoms->cFREEZE[i];
1170 for (m = 0; m < DIM; m++)
1172 if (!inputrec->opts.nFreeze[gf][m])
1174 pm[i][m] = sfm[i][m] + beta * pm[i][m];
1175 gpa -= pm[i][m] * sfm[i][m];
1176 /* f is negative gradient, thus the sign */
1185 /* Sum the gradient along the line across CPUs */
1188 gmx_sumd(1, &gpa, cr);
1191 /* Calculate the norm of the search vector */
1192 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1194 /* Just in case stepsize reaches zero due to numerical precision... */
1197 stepsize = inputrec->em_stepsize / pnorm;
1201 * Double check the value of the derivative in the search direction.
1202 * If it is positive it must be due to the old information in the
1203 * CG formula, so just remove that and start over with beta=0.
1204 * This corresponds to a steepest descent step.
1209 step--; /* Don't count this step since we are restarting */
1210 continue; /* Go back to the beginning of the big for-loop */
1213 /* Calculate minimum allowed stepsize, before the average (norm)
1214 * relative change in coordinate is smaller than precision
1217 auto s_min_x = makeArrayRef(s_min->s.x);
1218 for (int i = 0; i < mdatoms->homenr; i++)
1220 for (m = 0; m < DIM; m++)
1222 tmp = fabs(s_min_x[i][m]);
1227 tmp = pm[i][m] / tmp;
1228 minstep += tmp * tmp;
1231 /* Add up from all CPUs */
1234 gmx_sumd(1, &minstep, cr);
1237 minstep = GMX_REAL_EPS / sqrt(minstep / (3 * top_global->natoms));
1239 if (stepsize < minstep)
1245 /* Write coordinates if necessary */
1246 do_x = do_per_step(step, inputrec->nstxout);
1247 do_f = do_per_step(step, inputrec->nstfout);
1249 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, step, s_min,
1250 state_global, observablesHistory);
1252 /* Take a step downhill.
1253 * In theory, we should minimize the function along this direction.
1254 * That is quite possible, but it turns out to take 5-10 function evaluations
1255 * for each line. However, we dont really need to find the exact minimum -
1256 * it is much better to start a new CG step in a modified direction as soon
1257 * as we are close to it. This will save a lot of energy evaluations.
1259 * In practice, we just try to take a single step.
1260 * If it worked (i.e. lowered the energy), we increase the stepsize but
1261 * the continue straight to the next CG step without trying to find any minimum.
1262 * If it didn't work (higher energy), there must be a minimum somewhere between
1263 * the old position and the new one.
1265 * Due to the finite numerical accuracy, it turns out that it is a good idea
1266 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1267 * This leads to lower final energies in the tests I've done. / Erik
1269 s_a->epot = s_min->epot;
1271 c = a + stepsize; /* reference position along line is zero */
1273 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1275 em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec, imdSession,
1276 pull_work, s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
1279 /* Take a trial step (new coords in s_c) */
1280 do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c, constr, -1);
1283 /* Calculate energy for the trial step */
1284 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1286 /* Calc derivative along line */
1287 const rvec* pc = s_c->s.cg_p.rvec_array();
1288 const rvec* sfc = s_c->f.rvec_array();
1290 for (int i = 0; i < mdatoms->homenr; i++)
1292 for (m = 0; m < DIM; m++)
1294 gpc -= pc[i][m] * sfc[i][m]; /* f is negative gradient, thus the sign */
1297 /* Sum the gradient along the line across CPUs */
1300 gmx_sumd(1, &gpc, cr);
1303 /* This is the max amount of increase in energy we tolerate */
1304 tmp = std::sqrt(GMX_REAL_EPS) * fabs(s_a->epot);
1306 /* Accept the step if the energy is lower, or if it is not significantly higher
1307 * and the line derivative is still negative.
1309 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1312 /* Great, we found a better energy. Increase step for next iteration
1313 * if we are still going down, decrease it otherwise
1317 stepsize *= 1.618034; /* The golden section */
1321 stepsize *= 0.618034; /* 1/golden section */
1326 /* New energy is the same or higher. We will have to do some work
1327 * to find a smaller value in the interval. Take smaller step next time!
1330 stepsize *= 0.618034;
1334 /* OK, if we didn't find a lower value we will have to locate one now - there must
1335 * be one in the interval [a=0,c].
1336 * The same thing is valid here, though: Don't spend dozens of iterations to find
1337 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1338 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1340 * I also have a safeguard for potentially really pathological functions so we never
1341 * take more than 20 steps before we give up ...
1343 * If we already found a lower value we just skip this step and continue to the update.
1352 /* Select a new trial point.
1353 * If the derivatives at points a & c have different sign we interpolate to zero,
1354 * otherwise just do a bisection.
1356 if (gpa < 0 && gpc > 0)
1358 b = a + gpa * (a - c) / (gpc - gpa);
1365 /* safeguard if interpolation close to machine accuracy causes errors:
1366 * never go outside the interval
1368 if (b <= a || b >= c)
1373 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1375 /* Reload the old state */
1376 em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec, imdSession, pull_work,
1377 s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
1380 /* Take a trial step to this new point - new coords in s_b */
1381 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b, constr, -1);
1384 /* Calculate energy for the trial step */
1385 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1387 /* p does not change within a step, but since the domain decomposition
1388 * might change, we have to use cg_p of s_b here.
1390 const rvec* pb = s_b->s.cg_p.rvec_array();
1391 const rvec* sfb = s_b->f.rvec_array();
1393 for (int i = 0; i < mdatoms->homenr; i++)
1395 for (m = 0; m < DIM; m++)
1397 gpb -= pb[i][m] * sfb[i][m]; /* f is negative gradient, thus the sign */
1400 /* Sum the gradient along the line across CPUs */
1403 gmx_sumd(1, &gpb, cr);
1408 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a->epot, s_b->epot,
1412 epot_repl = s_b->epot;
1414 /* Keep one of the intervals based on the value of the derivative at the new point */
1417 /* Replace c endpoint with b */
1418 swap_em_state(&s_b, &s_c);
1424 /* Replace a endpoint with b */
1425 swap_em_state(&s_b, &s_a);
1431 * Stop search as soon as we find a value smaller than the endpoints.
1432 * Never run more than 20 steps, no matter what.
1435 } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) && (nminstep < 20));
1437 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot) * GMX_REAL_EPS || nminstep >= 20)
1439 /* OK. We couldn't find a significantly lower energy.
1440 * If beta==0 this was steepest descent, and then we give up.
1441 * If not, set beta=0 and restart with steepest descent before quitting.
1451 /* Reset memory before giving up */
1457 /* Select min energy state of A & C, put the best in B.
1459 if (s_c->epot < s_a->epot)
1463 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c->epot,
1466 swap_em_state(&s_b, &s_c);
1473 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a->epot,
1476 swap_em_state(&s_b, &s_a);
1484 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n", s_c->epot);
1486 swap_em_state(&s_b, &s_c);
1490 /* new search direction */
1491 /* beta = 0 means forget all memory and restart with steepest descents. */
1492 if (nstcg && ((step % nstcg) == 0))
1498 /* s_min->fnorm cannot be zero, because then we would have converged
1502 /* Polak-Ribiere update.
1503 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1505 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1507 /* Limit beta to prevent oscillations */
1508 if (fabs(beta) > 5.0)
1514 /* update positions */
1515 swap_em_state(&s_min, &s_b);
1518 /* Print it if necessary */
1521 if (mdrunOptions.verbose)
1523 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1524 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step,
1525 s_min->epot, s_min->fnorm / sqrtNumAtoms, s_min->fmax, s_min->a_fmax + 1);
1528 /* Store the new (lower) energies */
1529 matrix nullBox = {};
1530 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1531 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1532 nullptr, vir, pres, nullptr, mu_tot, constr);
1534 do_log = do_per_step(step, inputrec->nstlog);
1535 do_ene = do_per_step(step, inputrec->nstenergy);
1537 imdSession->fillEnergyRecord(step, TRUE);
1541 energyOutput.printHeader(fplog, step, step);
1543 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1544 do_log ? fplog : nullptr, step, step, fcd, nullptr);
1547 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1548 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0))
1550 imdSession->sendPositionsAndEnergies();
1553 /* Stop when the maximum force lies below tolerance.
1554 * If we have reached machine precision, converged is already set to true.
1556 converged = converged || (s_min->fmax < inputrec->em_tol);
1558 } /* End of the loop */
1562 step--; /* we never took that last step in this case */
1564 if (s_min->fmax > inputrec->em_tol)
1568 warn_step(fplog, inputrec->em_tol, s_min->fmax, step - 1 == number_steps, FALSE);
1575 /* If we printed energy and/or logfile last step (which was the last step)
1576 * we don't have to do it again, but otherwise print the final values.
1580 /* Write final value to log since we didn't do anything the last step */
1581 energyOutput.printHeader(fplog, step, step);
1583 if (!do_ene || !do_log)
1585 /* Write final energy file entries */
1586 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1587 !do_log ? fplog : nullptr, step, step, fcd, nullptr);
1591 /* Print some stuff... */
1594 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1598 * For accurate normal mode calculation it is imperative that we
1599 * store the last conformation into the full precision binary trajectory.
1601 * However, we should only do it if we did NOT already write this step
1602 * above (which we did if do_x or do_f was true).
1604 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1605 * in the trajectory with the same step number.
1607 do_x = !do_per_step(step, inputrec->nstxout);
1608 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1610 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec,
1611 step, s_min, state_global, observablesHistory);
1616 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1617 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1618 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1620 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1623 finish_em(cr, outf, walltime_accounting, wcycle);
1625 /* To print the actual number of steps we needed somewhere */
1626 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1630 void LegacySimulator::do_lbfgs()
1632 static const char* LBFGS = "Low-Memory BFGS Minimizer";
1635 gmx_global_stat_t gstat;
1637 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1638 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1639 real * rho, *alpha, *p, *s, **dx, **dg;
1640 real a, b, c, maxdelta, delta;
1642 real dgdx, dgdg, sq, yr, beta;
1644 rvec mu_tot = { 0 };
1645 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1647 int start, end, number_steps;
1648 int i, k, m, n, gf, step;
1650 auto mdatoms = mdAtoms->mdatoms();
1655 "Note that activating L-BFGS energy minimization via the "
1656 "integrator .mdp option and the command gmx mdrun may "
1657 "be available in a different form in a future version of GROMACS, "
1658 "e.g. gmx minimize and an .mdp option.");
1662 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1665 if (nullptr != constr)
1669 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1670 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1673 n = 3 * state_global->natoms;
1674 nmaxcorr = inputrec->nbfgscorr;
1679 snew(rho, nmaxcorr);
1680 snew(alpha, nmaxcorr);
1683 for (i = 0; i < nmaxcorr; i++)
1689 for (i = 0; i < nmaxcorr; i++)
1698 init_em(fplog, mdlog, LBFGS, cr, inputrec, imdSession, pull_work, state_global, top_global,
1699 &ems, &top, nrnb, fr, &graph, mdAtoms, &gstat, vsite, constr, nullptr);
1701 init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier,
1702 inputrec, top_global, nullptr, wcycle, StartingBehavior::NewSimulation);
1703 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
1704 false, StartingBehavior::NewSimulation, mdModulesNotifier);
1707 end = mdatoms->homenr;
1709 /* We need 4 working states */
1710 em_state_t s0{}, s1{}, s2{}, s3{};
1711 em_state_t* sa = &s0;
1712 em_state_t* sb = &s1;
1713 em_state_t* sc = &s2;
1714 em_state_t* last = &s3;
1715 /* Initialize by copying the state from ems (we could skip x and f here) */
1720 /* Print to log file */
1721 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1723 do_log = do_ene = do_x = do_f = TRUE;
1725 /* Max number of steps */
1726 number_steps = inputrec->nsteps;
1728 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1730 for (i = start; i < end; i++)
1732 if (mdatoms->cFREEZE)
1734 gf = mdatoms->cFREEZE[i];
1736 for (m = 0; m < DIM; m++)
1738 frozen[3 * i + m] = (inputrec->opts.nFreeze[gf][m] != 0);
1743 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1747 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1752 construct_vsites(vsite, state_global->x.rvec_array(), 1, nullptr, top.idef.iparams,
1753 top.idef.il, fr->ePBC, fr->bMolPBC, cr, state_global->box);
1756 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1757 /* do_force always puts the charge groups in the box and shifts again
1758 * We do not unshift, so molecules are always whole
1761 EnergyEvaluator energyEvaluator{
1762 fplog, mdlog, cr, ms, top_global, &top, inputrec,
1763 imdSession, pull_work, nrnb, wcycle, gstat, vsite, constr,
1764 fcd, graph, mdAtoms, fr, runScheduleWork, enerd
1766 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1770 /* Copy stuff to the energy bin for easy printing etc. */
1771 matrix nullBox = {};
1772 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
1773 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
1774 nullptr, vir, pres, nullptr, mu_tot, constr);
1776 energyOutput.printHeader(fplog, step, step);
1777 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step,
1778 step, fcd, nullptr);
1781 /* Set the initial step.
1782 * since it will be multiplied by the non-normalized search direction
1783 * vector (force vector the first time), we scale it by the
1784 * norm of the force.
1789 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1790 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1791 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1792 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
1793 fprintf(stderr, "\n");
1794 /* and copy to the log file too... */
1795 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1796 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1797 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
1798 fprintf(fplog, "\n");
1801 // Point is an index to the memory of search directions, where 0 is the first one.
1804 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1805 real* fInit = static_cast<real*>(ems.f.rvec_array()[0]);
1806 for (i = 0; i < n; i++)
1810 dx[point][i] = fInit[i]; /* Initial search direction */
1818 // Stepsize will be modified during the search, and actually it is not critical
1819 // (the main efficiency in the algorithm comes from changing directions), but
1820 // we still need an initial value, so estimate it as the inverse of the norm
1821 // so we take small steps where the potential fluctuates a lot.
1822 stepsize = 1.0 / ems.fnorm;
1824 /* Start the loop over BFGS steps.
1825 * Each successful step is counted, and we continue until
1826 * we either converge or reach the max number of steps.
1831 /* Set the gradient from the force */
1833 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1836 /* Write coordinates if necessary */
1837 do_x = do_per_step(step, inputrec->nstxout);
1838 do_f = do_per_step(step, inputrec->nstfout);
1843 mdof_flags |= MDOF_X;
1848 mdof_flags |= MDOF_F;
1853 mdof_flags |= MDOF_IMD;
1856 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global->natoms, step,
1857 static_cast<real>(step), &ems.s, state_global,
1858 observablesHistory, ems.f);
1860 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1862 /* make s a pointer to current search direction - point=0 first time we get here */
1865 real* xx = static_cast<real*>(ems.s.x.rvec_array()[0]);
1866 real* ff = static_cast<real*>(ems.f.rvec_array()[0]);
1868 // calculate line gradient in position A
1869 for (gpa = 0, i = 0; i < n; i++)
1871 gpa -= s[i] * ff[i];
1874 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1875 * relative change in coordinate is smaller than precision
1877 for (minstep = 0, i = 0; i < n; i++)
1885 minstep += tmp * tmp;
1887 minstep = GMX_REAL_EPS / sqrt(minstep / n);
1889 if (stepsize < minstep)
1895 // Before taking any steps along the line, store the old position
1897 real* lastx = static_cast<real*>(last->s.x.data()[0]);
1898 real* lastf = static_cast<real*>(last->f.data()[0]);
1903 /* Take a step downhill.
1904 * In theory, we should find the actual minimum of the function in this
1905 * direction, somewhere along the line.
1906 * That is quite possible, but it turns out to take 5-10 function evaluations
1907 * for each line. However, we dont really need to find the exact minimum -
1908 * it is much better to start a new BFGS step in a modified direction as soon
1909 * as we are close to it. This will save a lot of energy evaluations.
1911 * In practice, we just try to take a single step.
1912 * If it worked (i.e. lowered the energy), we increase the stepsize but
1913 * continue straight to the next BFGS step without trying to find any minimum,
1914 * i.e. we change the search direction too. If the line was smooth, it is
1915 * likely we are in a smooth region, and then it makes sense to take longer
1916 * steps in the modified search direction too.
1918 * If it didn't work (higher energy), there must be a minimum somewhere between
1919 * the old position and the new one. Then we need to start by finding a lower
1920 * value before we change search direction. Since the energy was apparently
1921 * quite rough, we need to decrease the step size.
1923 * Due to the finite numerical accuracy, it turns out that it is a good idea
1924 * to accept a SMALL increase in energy, if the derivative is still downhill.
1925 * This leads to lower final energies in the tests I've done. / Erik
1928 // State "A" is the first position along the line.
1929 // reference position along line is initially zero
1932 // Check stepsize first. We do not allow displacements
1933 // larger than emstep.
1937 // Pick a new position C by adding stepsize to A.
1940 // Calculate what the largest change in any individual coordinate
1941 // would be (translation along line * gradient along line)
1943 for (i = 0; i < n; i++)
1946 if (delta > maxdelta)
1951 // If any displacement is larger than the stepsize limit, reduce the step
1952 if (maxdelta > inputrec->em_stepsize)
1956 } while (maxdelta > inputrec->em_stepsize);
1958 // Take a trial step and move the coordinate array xc[] to position C
1959 real* xc = static_cast<real*>(sc->s.x.rvec_array()[0]);
1960 for (i = 0; i < n; i++)
1962 xc[i] = lastx[i] + c * s[i];
1966 // Calculate energy for the trial step in position C
1967 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
1969 // Calc line gradient in position C
1970 real* fc = static_cast<real*>(sc->f.rvec_array()[0]);
1971 for (gpc = 0, i = 0; i < n; i++)
1973 gpc -= s[i] * fc[i]; /* f is negative gradient, thus the sign */
1975 /* Sum the gradient along the line across CPUs */
1978 gmx_sumd(1, &gpc, cr);
1981 // This is the max amount of increase in energy we tolerate.
1982 // By allowing VERY small changes (close to numerical precision) we
1983 // frequently find even better (lower) final energies.
1984 tmp = std::sqrt(GMX_REAL_EPS) * fabs(sa->epot);
1986 // Accept the step if the energy is lower in the new position C (compared to A),
1987 // or if it is not significantly higher and the line derivative is still negative.
1988 foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
1989 // If true, great, we found a better energy. We no longer try to alter the
1990 // stepsize, but simply accept this new better position. The we select a new
1991 // search direction instead, which will be much more efficient than continuing
1992 // to take smaller steps along a line. Set fnorm based on the new C position,
1993 // which will be used to update the stepsize to 1/fnorm further down.
1995 // If false, the energy is NOT lower in point C, i.e. it will be the same
1996 // or higher than in point A. In this case it is pointless to move to point C,
1997 // so we will have to do more iterations along the same line to find a smaller
1998 // value in the interval [A=0.0,C].
1999 // Here, A is still 0.0, but that will change when we do a search in the interval
2000 // [0.0,C] below. That search we will do by interpolation or bisection rather
2001 // than with the stepsize, so no need to modify it. For the next search direction
2002 // it will be reset to 1/fnorm anyway.
2006 // OK, if we didn't find a lower value we will have to locate one now - there must
2007 // be one in the interval [a,c].
2008 // The same thing is valid here, though: Don't spend dozens of iterations to find
2009 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2010 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2011 // I also have a safeguard for potentially really pathological functions so we never
2012 // take more than 20 steps before we give up.
2013 // If we already found a lower value we just skip this step and continue to the update.
2018 // Select a new trial point B in the interval [A,C].
2019 // If the derivatives at points a & c have different sign we interpolate to zero,
2020 // otherwise just do a bisection since there might be multiple minima/maxima
2021 // inside the interval.
2022 if (gpa < 0 && gpc > 0)
2024 b = a + gpa * (a - c) / (gpc - gpa);
2031 /* safeguard if interpolation close to machine accuracy causes errors:
2032 * never go outside the interval
2034 if (b <= a || b >= c)
2039 // Take a trial step to point B
2040 real* xb = static_cast<real*>(sb->s.x.rvec_array()[0]);
2041 for (i = 0; i < n; i++)
2043 xb[i] = lastx[i] + b * s[i];
2047 // Calculate energy for the trial step in point B
2048 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2051 // Calculate gradient in point B
2052 real* fb = static_cast<real*>(sb->f.rvec_array()[0]);
2053 for (gpb = 0, i = 0; i < n; i++)
2055 gpb -= s[i] * fb[i]; /* f is negative gradient, thus the sign */
2057 /* Sum the gradient along the line across CPUs */
2060 gmx_sumd(1, &gpb, cr);
2063 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2064 // at the new point B, and rename the endpoints of this new interval A and C.
2067 /* Replace c endpoint with b */
2069 /* copy state b to c */
2074 /* Replace a endpoint with b */
2076 /* copy state b to a */
2081 * Stop search as soon as we find a value smaller than the endpoints,
2082 * or if the tolerance is below machine precision.
2083 * Never run more than 20 steps, no matter what.
2086 } while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2088 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2090 /* OK. We couldn't find a significantly lower energy.
2091 * If ncorr==0 this was steepest descent, and then we give up.
2092 * If not, reset memory to restart as steepest descent before quitting.
2104 /* Search in gradient direction */
2105 for (i = 0; i < n; i++)
2107 dx[point][i] = ff[i];
2109 /* Reset stepsize */
2110 stepsize = 1.0 / fnorm;
2115 /* Select min energy state of A & C, put the best in xx/ff/Epot
2117 if (sc->epot < sa->epot)
2138 /* Update the memory information, and calculate a new
2139 * approximation of the inverse hessian
2142 /* Have new data in Epot, xx, ff */
2143 if (ncorr < nmaxcorr)
2148 for (i = 0; i < n; i++)
2150 dg[point][i] = lastf[i] - ff[i];
2151 dx[point][i] *= step_taken;
2156 for (i = 0; i < n; i++)
2158 dgdg += dg[point][i] * dg[point][i];
2159 dgdx += dg[point][i] * dx[point][i];
2164 rho[point] = 1.0 / dgdx;
2167 if (point >= nmaxcorr)
2173 for (i = 0; i < n; i++)
2180 /* Recursive update. First go back over the memory points */
2181 for (k = 0; k < ncorr; k++)
2190 for (i = 0; i < n; i++)
2192 sq += dx[cp][i] * p[i];
2195 alpha[cp] = rho[cp] * sq;
2197 for (i = 0; i < n; i++)
2199 p[i] -= alpha[cp] * dg[cp][i];
2203 for (i = 0; i < n; i++)
2208 /* And then go forward again */
2209 for (k = 0; k < ncorr; k++)
2212 for (i = 0; i < n; i++)
2214 yr += p[i] * dg[cp][i];
2217 beta = rho[cp] * yr;
2218 beta = alpha[cp] - beta;
2220 for (i = 0; i < n; i++)
2222 p[i] += beta * dx[cp][i];
2232 for (i = 0; i < n; i++)
2236 dx[point][i] = p[i];
2244 /* Print it if necessary */
2247 if (mdrunOptions.verbose)
2249 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2250 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step,
2251 ems.epot, ems.fnorm / sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2254 /* Store the new (lower) energies */
2255 matrix nullBox = {};
2256 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step), mdatoms->tmass,
2257 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
2258 nullptr, vir, pres, nullptr, mu_tot, constr);
2260 do_log = do_per_step(step, inputrec->nstlog);
2261 do_ene = do_per_step(step, inputrec->nstenergy);
2263 imdSession->fillEnergyRecord(step, TRUE);
2267 energyOutput.printHeader(fplog, step, step);
2269 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2270 do_log ? fplog : nullptr, step, step, fcd, nullptr);
2273 /* Send x and E to IMD client, if bIMD is TRUE. */
2274 if (imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0) && MASTER(cr))
2276 imdSession->sendPositionsAndEnergies();
2279 // Reset stepsize in we are doing more iterations
2282 /* Stop when the maximum force lies below tolerance.
2283 * If we have reached machine precision, converged is already set to true.
2285 converged = converged || (ems.fmax < inputrec->em_tol);
2287 } /* End of the loop */
2291 step--; /* we never took that last step in this case */
2293 if (ems.fmax > inputrec->em_tol)
2297 warn_step(fplog, inputrec->em_tol, ems.fmax, step - 1 == number_steps, FALSE);
2302 /* If we printed energy and/or logfile last step (which was the last step)
2303 * we don't have to do it again, but otherwise print the final values.
2305 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2307 energyOutput.printHeader(fplog, step, step);
2309 if (!do_ene || !do_log) /* Write final energy file entries */
2311 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2312 !do_log ? fplog : nullptr, step, step, fcd, nullptr);
2315 /* Print some stuff... */
2318 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2322 * For accurate normal mode calculation it is imperative that we
2323 * store the last conformation into the full precision binary trajectory.
2325 * However, we should only do it if we did NOT already write this step
2326 * above (which we did if do_x or do_f was true).
2328 do_x = !do_per_step(step, inputrec->nstxout);
2329 do_f = !do_per_step(step, inputrec->nstfout);
2330 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec,
2331 step, &ems, state_global, observablesHistory);
2335 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2336 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2337 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2339 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2342 finish_em(cr, outf, walltime_accounting, wcycle);
2344 /* To print the actual number of steps we needed somewhere */
2345 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2348 void LegacySimulator::do_steep()
2350 const char* SD = "Steepest Descents";
2352 gmx_global_stat_t gstat;
2356 gmx_bool bDone, bAbort, do_x, do_f;
2358 rvec mu_tot = { 0 };
2361 int steps_accepted = 0;
2362 auto mdatoms = mdAtoms->mdatoms();
2367 "Note that activating steepest-descent energy minimization via the "
2368 "integrator .mdp option and the command gmx mdrun may "
2369 "be available in a different form in a future version of GROMACS, "
2370 "e.g. gmx minimize and an .mdp option.");
2372 /* Create 2 states on the stack and extract pointers that we will swap */
2373 em_state_t s0{}, s1{};
2374 em_state_t* s_min = &s0;
2375 em_state_t* s_try = &s1;
2377 /* Init em and store the local state in s_try */
2378 init_em(fplog, mdlog, SD, cr, inputrec, imdSession, pull_work, state_global, top_global, s_try,
2379 &top, nrnb, fr, &graph, mdAtoms, &gstat, vsite, constr, nullptr);
2381 init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier,
2382 inputrec, top_global, nullptr, wcycle, StartingBehavior::NewSimulation);
2383 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr,
2384 false, StartingBehavior::NewSimulation, mdModulesNotifier);
2386 /* Print to log file */
2387 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2389 /* Set variables for stepsize (in nm). This is the largest
2390 * step that we are going to make in any direction.
2392 ustep = inputrec->em_stepsize;
2395 /* Max number of steps */
2396 nsteps = inputrec->nsteps;
2400 /* Print to the screen */
2401 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2405 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2407 EnergyEvaluator energyEvaluator{
2408 fplog, mdlog, cr, ms, top_global, &top, inputrec,
2409 imdSession, pull_work, nrnb, wcycle, gstat, vsite, constr,
2410 fcd, graph, mdAtoms, fr, runScheduleWork, enerd
2413 /**** HERE STARTS THE LOOP ****
2414 * count is the counter for the number of steps
2415 * bDone will be TRUE when the minimization has converged
2416 * bAbort will be TRUE when nsteps steps have been performed or when
2417 * the stepsize becomes smaller than is reasonable for machine precision
2422 while (!bDone && !bAbort)
2424 bAbort = (nsteps >= 0) && (count == nsteps);
2426 /* set new coordinates, except for first step */
2427 bool validStep = true;
2430 validStep = do_em_step(cr, inputrec, mdatoms, s_min, stepsize, &s_min->f, s_try, constr, count);
2435 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2439 // Signal constraint error during stepping with energy=inf
2440 s_try->epot = std::numeric_limits<real>::infinity();
2445 energyOutput.printHeader(fplog, count, count);
2450 s_min->epot = s_try->epot;
2453 /* Print it if necessary */
2456 if (mdrunOptions.verbose)
2458 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2459 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax + 1,
2460 ((count == 0) || (s_try->epot < s_min->epot)) ? '\n' : '\r');
2464 if ((count == 0) || (s_try->epot < s_min->epot))
2466 /* Store the new (lower) energies */
2467 matrix nullBox = {};
2468 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(count), mdatoms->tmass,
2469 enerd, nullptr, nullptr, nullptr, nullBox, nullptr,
2470 nullptr, vir, pres, nullptr, mu_tot, constr);
2472 imdSession->fillEnergyRecord(count, TRUE);
2474 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2475 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2476 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, do_dr, do_or,
2477 fplog, count, count, fcd, nullptr);
2482 /* Now if the new energy is smaller than the previous...
2483 * or if this is the first step!
2484 * or if we did random steps!
2487 if ((count == 0) || (s_try->epot < s_min->epot))
2491 /* Test whether the convergence criterion is met... */
2492 bDone = (s_try->fmax < inputrec->em_tol);
2494 /* Copy the arrays for force, positions and energy */
2495 /* The 'Min' array always holds the coords and forces of the minimal
2497 swap_em_state(&s_min, &s_try);
2503 /* Write to trn, if necessary */
2504 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2505 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2506 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, count, s_min,
2507 state_global, observablesHistory);
2511 /* If energy is not smaller make the step smaller... */
2514 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2516 /* Reload the old state */
2517 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession,
2518 pull_work, s_min, &top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
2522 // If the force is very small after finishing minimization,
2523 // we risk dividing by zero when calculating the step size.
2524 // So we check first if the minimization has stopped before
2525 // trying to obtain a new step size.
2528 /* Determine new step */
2529 stepsize = ustep / s_min->fmax;
2532 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2534 if (count == nsteps || ustep < 1e-12)
2536 if (count == nsteps || ustep < 1e-6)
2541 warn_step(fplog, inputrec->em_tol, s_min->fmax, count == nsteps, constr != nullptr);
2546 /* Send IMD energies and positions, if bIMD is TRUE. */
2547 if (imdSession->run(count, TRUE, state_global->box,
2548 MASTER(cr) ? state_global->x.rvec_array() : nullptr, 0)
2551 imdSession->sendPositionsAndEnergies();
2555 } /* End of the loop */
2557 /* Print some data... */
2560 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2562 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2563 top_global, inputrec, count, s_min, state_global, observablesHistory);
2567 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2569 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2570 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2573 finish_em(cr, outf, walltime_accounting, wcycle);
2575 /* To print the actual number of steps we needed somewhere */
2576 inputrec->nsteps = count;
2578 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2581 void LegacySimulator::do_nm()
2583 const char* NM = "Normal Mode Analysis";
2586 gmx_global_stat_t gstat;
2589 rvec mu_tot = { 0 };
2591 gmx_bool bSparse; /* use sparse matrix storage format */
2593 gmx_sparsematrix_t* sparse_matrix = nullptr;
2594 real* full_matrix = nullptr;
2596 /* added with respect to mdrun */
2598 real der_range = 10.0 * std::sqrt(GMX_REAL_EPS);
2600 bool bIsMaster = MASTER(cr);
2601 auto mdatoms = mdAtoms->mdatoms();
2606 "Note that activating normal-mode analysis via the integrator "
2607 ".mdp option and the command gmx mdrun may "
2608 "be available in a different form in a future version of GROMACS, "
2609 "e.g. gmx normal-modes.");
2611 if (constr != nullptr)
2615 "Constraints present with Normal Mode Analysis, this combination is not supported");
2618 gmx_shellfc_t* shellfc;
2620 em_state_t state_work{};
2622 /* Init em and store the local state in state_minimum */
2623 init_em(fplog, mdlog, NM, cr, inputrec, imdSession, pull_work, state_global, top_global,
2624 &state_work, &top, nrnb, fr, &graph, mdAtoms, &gstat, vsite, constr, &shellfc);
2626 init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier,
2627 inputrec, top_global, nullptr, wcycle, StartingBehavior::NewSimulation);
2629 std::vector<int> atom_index = get_atom_index(top_global);
2630 std::vector<gmx::RVec> fneg(atom_index.size(), { 0, 0, 0 });
2631 snew(dfdx, atom_index.size());
2637 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2638 " which MIGHT not be accurate enough for normal mode analysis.\n"
2639 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2640 " are fairly modest even if you recompile in double precision.\n\n");
2644 /* Check if we can/should use sparse storage format.
2646 * Sparse format is only useful when the Hessian itself is sparse, which it
2647 * will be when we use a cutoff.
2648 * For small systems (n<1000) it is easier to always use full matrix format, though.
2650 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2652 GMX_LOG(mdlog.warning)
2653 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2656 else if (atom_index.size() < 1000)
2658 GMX_LOG(mdlog.warning)
2659 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2665 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2669 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2670 sz = DIM * atom_index.size();
2672 fprintf(stderr, "Allocating Hessian memory...\n\n");
2676 sparse_matrix = gmx_sparsematrix_init(sz);
2677 sparse_matrix->compressed_symmetric = TRUE;
2681 snew(full_matrix, sz * sz);
2684 /* Write start time and temperature */
2685 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2687 /* fudge nr of steps to nr of atoms */
2688 inputrec->nsteps = atom_index.size() * 2;
2692 fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2693 *(top_global->name), inputrec->nsteps);
2696 nnodes = cr->nnodes;
2698 /* Make evaluate_energy do a single node force calculation */
2700 EnergyEvaluator energyEvaluator{
2701 fplog, mdlog, cr, ms, top_global, &top, inputrec,
2702 imdSession, pull_work, nrnb, wcycle, gstat, vsite, constr,
2703 fcd, graph, mdAtoms, fr, runScheduleWork, enerd
2705 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2706 cr->nnodes = nnodes;
2708 /* if forces are not small, warn user */
2709 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2711 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2712 if (state_work.fmax > 1.0e-3)
2714 GMX_LOG(mdlog.warning)
2716 "The force is probably not small enough to "
2717 "ensure that you are at a minimum.\n"
2718 "Be aware that negative eigenvalues may occur\n"
2719 "when the resulting matrix is diagonalized.");
2722 /***********************************************************
2724 * Loop over all pairs in matrix
2726 * do_force called twice. Once with positive and
2727 * once with negative displacement
2729 ************************************************************/
2731 /* Steps are divided one by one over the nodes */
2733 auto state_work_x = makeArrayRef(state_work.s.x);
2734 auto state_work_f = makeArrayRef(state_work.f);
2735 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
2737 size_t atom = atom_index[aid];
2738 for (size_t d = 0; d < DIM; d++)
2741 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2744 x_min = state_work_x[atom][d];
2746 for (unsigned int dx = 0; (dx < 2); dx++)
2750 state_work_x[atom][d] = x_min - der_range;
2754 state_work_x[atom][d] = x_min + der_range;
2757 /* Make evaluate_energy do a single node force calculation */
2761 /* Now is the time to relax the shells */
2762 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose, nullptr, step, inputrec,
2763 imdSession, pull_work, bNS, force_flags, &top, constr, enerd,
2764 fcd, state_work.s.natoms, state_work.s.x.arrayRefWithPadding(),
2765 state_work.s.v.arrayRefWithPadding(), state_work.s.box,
2766 state_work.s.lambda, &state_work.s.hist,
2767 state_work.f.arrayRefWithPadding(), vir, mdatoms, nrnb,
2768 wcycle, graph, shellfc, fr, runScheduleWork, t, mu_tot,
2769 vsite, DDBalanceRegionHandler(nullptr));
2775 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid * 2 + dx, FALSE);
2778 cr->nnodes = nnodes;
2782 std::copy(state_work_f.begin(), state_work_f.begin() + atom_index.size(),
2787 /* x is restored to original */
2788 state_work_x[atom][d] = x_min;
2790 for (size_t j = 0; j < atom_index.size(); j++)
2792 for (size_t k = 0; (k < DIM); k++)
2794 dfdx[j][k] = -(state_work_f[atom_index[j]][k] - fneg[j][k]) / (2 * der_range);
2801 # define mpi_type GMX_MPI_REAL
2802 MPI_Send(dfdx[0], atom_index.size() * DIM, mpi_type, MASTER(cr), cr->nodeid,
2803 cr->mpi_comm_mygroup);
2808 for (index node = 0; (node < nnodes && aid + node < ssize(atom_index)); node++)
2814 MPI_Recv(dfdx[0], atom_index.size() * DIM, mpi_type, node, node,
2815 cr->mpi_comm_mygroup, &stat);
2820 row = (aid + node) * DIM + d;
2822 for (size_t j = 0; j < atom_index.size(); j++)
2824 for (size_t k = 0; k < DIM; k++)
2830 if (col >= row && dfdx[j][k] != 0.0)
2832 gmx_sparsematrix_increment_value(sparse_matrix, row, col, dfdx[j][k]);
2837 full_matrix[row * sz + col] = dfdx[j][k];
2844 if (mdrunOptions.verbose && fplog)
2849 /* write progress */
2850 if (bIsMaster && mdrunOptions.verbose)
2852 fprintf(stderr, "\rFinished step %d out of %td",
2853 std::min<int>(atom + nnodes, atom_index.size()), ssize(atom_index));
2860 fprintf(stderr, "\n\nWriting Hessian...\n");
2861 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2864 finish_em(cr, outf, walltime_accounting, wcycle);
2866 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size() * 2);