2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 * \brief This file defines integrators for energy minimization
42 * \author Berk Hess <hess@kth.se>
43 * \author Erik Lindahl <erik@kth.se>
44 * \ingroup module_mdrun
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/collect.h"
59 #include "gromacs/domdec/dlbtiming.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/mdsetup.h"
63 #include "gromacs/domdec/partition.h"
64 #include "gromacs/ewald/pme_pp.h"
65 #include "gromacs/fileio/confio.h"
66 #include "gromacs/fileio/mtxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/imd/imd.h"
70 #include "gromacs/linearalgebra/sparsematrix.h"
71 #include "gromacs/listed_forces/listed_forces.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/mdlib/constr.h"
75 #include "gromacs/mdlib/coupling.h"
76 #include "gromacs/mdlib/dispersioncorrection.h"
77 #include "gromacs/mdlib/ebin.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/energyoutput.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/gmx_omp_nthreads.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/stat.h"
87 #include "gromacs/mdlib/tgroup.h"
88 #include "gromacs/mdlib/trajectory_writing.h"
89 #include "gromacs/mdlib/update.h"
90 #include "gromacs/mdlib/vsite.h"
91 #include "gromacs/mdrunutility/handlerestart.h"
92 #include "gromacs/mdrunutility/printtime.h"
93 #include "gromacs/mdtypes/checkpointdata.h"
94 #include "gromacs/mdtypes/commrec.h"
95 #include "gromacs/mdtypes/forcebuffers.h"
96 #include "gromacs/mdtypes/forcerec.h"
97 #include "gromacs/mdtypes/inputrec.h"
98 #include "gromacs/mdtypes/interaction_const.h"
99 #include "gromacs/mdtypes/md_enums.h"
100 #include "gromacs/mdtypes/mdatom.h"
101 #include "gromacs/mdtypes/mdrunoptions.h"
102 #include "gromacs/mdtypes/state.h"
103 #include "gromacs/pbcutil/pbc.h"
104 #include "gromacs/timing/wallcycle.h"
105 #include "gromacs/timing/walltime_accounting.h"
106 #include "gromacs/topology/mtop_util.h"
107 #include "gromacs/topology/topology.h"
108 #include "gromacs/utility/cstringutil.h"
109 #include "gromacs/utility/exceptions.h"
110 #include "gromacs/utility/fatalerror.h"
111 #include "gromacs/utility/logger.h"
112 #include "gromacs/utility/smalloc.h"
114 #include "legacysimulator.h"
118 using gmx::MdrunScheduleWorkload;
120 using gmx::VirtualSitesHandler;
122 //! Utility structure for manipulating states during EM
123 typedef struct em_state
125 //! Copy of the global state
131 //! Norm of the force
139 //! Print the EM starting conditions
140 static void print_em_start(FILE* fplog,
142 gmx_walltime_accounting_t walltime_accounting,
143 gmx_wallcycle_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",
241 gmx_step_str(count, buf),
246 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n", alg, ftol, gmx_step_str(count, buf));
250 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
251 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
252 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm / sqrtNumAtoms);
254 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
255 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
256 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm / sqrtNumAtoms);
260 //! Compute the norm and max of the force array in parallel
261 static void get_f_norm_max(const t_commrec* cr,
262 const t_grpopts* opts,
264 gmx::ArrayRef<const gmx::RVec> f,
271 int la_max, a_max, start, end, i, m, gf;
273 /* This routine finds the largest force and returns it.
274 * On parallel machines the global max is taken.
280 end = mdatoms->homenr;
281 if (mdatoms->cFREEZE)
283 for (i = start; i < end; i++)
285 gf = mdatoms->cFREEZE[i];
287 for (m = 0; m < DIM; m++)
289 if (!opts->nFreeze[gf][m])
291 fam += gmx::square(f[i][m]);
304 for (i = start; i < end; i++)
316 if (la_max >= 0 && DOMAINDECOMP(cr))
318 a_max = cr->dd->globalAtomIndices[la_max];
326 snew(sum, 2 * cr->nnodes + 1);
327 sum[2 * cr->nodeid] = fmax2;
328 sum[2 * cr->nodeid + 1] = a_max;
329 sum[2 * cr->nnodes] = fnorm2;
330 gmx_sumd(2 * cr->nnodes + 1, sum, cr);
331 fnorm2 = sum[2 * cr->nnodes];
332 /* Determine the global maximum */
333 for (i = 0; i < cr->nnodes; i++)
335 if (sum[2 * i] > fmax2)
338 a_max = gmx::roundToInt(sum[2 * i + 1]);
346 *fnorm = sqrt(fnorm2);
358 //! Compute the norm of the force
359 static void get_state_f_norm_max(const t_commrec* cr, const t_grpopts* opts, t_mdatoms* mdatoms, em_state_t* ems)
361 get_f_norm_max(cr, opts, mdatoms, ems->f.view().force(), &ems->fnorm, &ems->fmax, &ems->a_fmax);
364 //! Initialize the energy minimization
365 static void init_em(FILE* fplog,
366 const gmx::MDLogger& mdlog,
369 const t_inputrec* ir,
370 gmx::ImdSession* imdSession,
372 t_state* state_global,
373 const gmx_mtop_t& top_global,
378 gmx::MDAtoms* mdAtoms,
379 gmx_global_stat_t* gstat,
380 VirtualSitesHandler* vsite,
381 gmx::Constraints* constr,
382 gmx_shellfc_t** shellfc)
388 fprintf(fplog, "Initiating %s\n", title);
393 state_global->ngtc = 0;
395 int* fep_state = MASTER(cr) ? &state_global->fep_state : nullptr;
396 gmx::ArrayRef<real> lambda = MASTER(cr) ? state_global->lambda : gmx::ArrayRef<real>();
398 fplog, *ir, gmx::arrayRefFromArray(ir->opts.ref_t, ir->opts.ngtc), MASTER(cr), fep_state, lambda);
400 if (ir->eI == IntegrationAlgorithm::NM)
402 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
404 *shellfc = init_shell_flexcon(stdout,
406 constr ? constr->numFlexibleConstraints() : 0,
409 thisRankHasDuty(cr, DUTY_PME));
413 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI),
414 "This else currently only handles energy minimizers, consider if your algorithm "
415 "needs shell/flexible-constraint support");
417 /* With energy minimization, shells and flexible constraints are
418 * automatically minimized when treated like normal DOFS.
420 if (shellfc != nullptr)
426 if (DOMAINDECOMP(cr))
428 dd_init_local_state(*cr->dd, state_global, &ems->s);
430 /* Distribute the charge groups over the nodes from the master node */
431 dd_partition_system(fplog,
452 dd_store_state(*cr->dd, &ems->s);
456 state_change_natoms(state_global, state_global->natoms);
457 /* Just copy the state */
458 ems->s = *state_global;
459 state_change_natoms(&ems->s, ems->s.natoms);
461 mdAlgorithmsSetupAtomData(
462 cr, *ir, top_global, top, fr, &ems->f, mdAtoms, constr, vsite, shellfc ? *shellfc : nullptr);
465 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[FreeEnergyPerturbationCouplingType::Mass]);
469 // TODO how should this cross-module support dependency be managed?
470 if (ir->eConstrAlg == ConstraintAlgorithm::Shake && gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
473 "Can not do energy minimization with %s, use %s\n",
474 enumValueToString(ConstraintAlgorithm::Shake),
475 enumValueToString(ConstraintAlgorithm::Lincs));
478 if (!ir->bContinuation)
480 /* Constrain the starting coordinates */
481 bool needsLogging = true;
482 bool computeEnergy = true;
483 bool computeVirial = false;
485 constr->apply(needsLogging,
490 ems->s.x.arrayRefWithPadding(),
491 ems->s.x.arrayRefWithPadding(),
494 ems->s.lambda[FreeEnergyPerturbationCouplingType::Fep],
496 gmx::ArrayRefWithPadding<RVec>(),
499 gmx::ConstraintVariable::Positions);
505 *gstat = global_stat_init(ir);
512 calc_shifts(ems->s.box, fr->shift_vec);
515 //! Finalize the minimization
516 static void finish_em(const t_commrec* cr,
518 gmx_walltime_accounting_t walltime_accounting,
519 gmx_wallcycle_t wcycle)
521 if (!thisRankHasDuty(cr, DUTY_PME))
523 /* Tell the PME only node to finish */
524 gmx_pme_send_finish(cr);
529 em_time_end(walltime_accounting, wcycle);
532 //! Swap two different EM states during minimization
533 static void swap_em_state(em_state_t** ems1, em_state_t** ems2)
542 //! Save the EM trajectory
543 static void write_em_traj(FILE* fplog,
549 const gmx_mtop_t& top_global,
550 const t_inputrec* ir,
553 t_state* state_global,
554 ObservablesHistory* observablesHistory)
560 mdof_flags |= MDOF_X;
564 mdof_flags |= MDOF_F;
567 /* If we want IMD output, set appropriate MDOF flag */
570 mdof_flags |= MDOF_IMD;
573 gmx::WriteCheckpointDataHolder checkpointDataHolder;
574 mdoutf_write_to_trajectory_files(fplog,
580 static_cast<double>(step),
584 state->f.view().force(),
585 &checkpointDataHolder);
587 if (confout != nullptr)
589 if (DOMAINDECOMP(cr))
591 /* If bX=true, x was collected to state_global in the call above */
594 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
596 cr->dd, state->s.ddp_count, state->s.ddp_count_cg_gl, state->s.cg_gl, state->s.x, globalXRef);
601 /* Copy the local state pointer */
602 state_global = &state->s;
607 if (ir->pbcType != PbcType::No && !ir->bPeriodicMols && DOMAINDECOMP(cr))
609 /* Make molecules whole only for confout writing */
610 do_pbc_mtop(ir->pbcType, state->s.box, &top_global, state_global->x.rvec_array());
613 write_sto_conf_mtop(confout,
616 state_global->x.rvec_array(),
624 //! \brief Do one minimization step
626 // \returns true when the step succeeded, false when a constraint error occurred
627 static bool do_em_step(const t_commrec* cr,
628 const t_inputrec* ir,
632 gmx::ArrayRefWithPadding<const gmx::RVec> force,
634 gmx::Constraints* constr,
641 int nthreads gmx_unused;
643 bool validStep = true;
648 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
650 gmx_incons("state mismatch in do_em_step");
653 s2->flags = s1->flags;
655 if (s2->natoms != s1->natoms)
657 state_change_natoms(s2, s1->natoms);
658 ems2->f.resize(s2->natoms);
660 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
662 s2->cg_gl.resize(s1->cg_gl.size());
665 copy_mat(s1->box, s2->box);
666 /* Copy free energy state */
667 s2->lambda = s1->lambda;
668 copy_mat(s1->box, s2->box);
673 nthreads = gmx_omp_nthreads_get(emntUpdate);
674 #pragma omp parallel num_threads(nthreads)
676 const rvec* x1 = s1->x.rvec_array();
677 rvec* x2 = s2->x.rvec_array();
678 const rvec* f = as_rvec_array(force.unpaddedArrayRef().data());
681 #pragma omp for schedule(static) nowait
682 for (int i = start; i < end; i++)
690 for (int m = 0; m < DIM; m++)
692 if (ir->opts.nFreeze[gf][m])
698 x2[i][m] = x1[i][m] + a * f[i][m];
702 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
705 if (s2->flags & enumValueToBitMask(StateEntry::Cgp))
707 /* Copy the CG p vector */
708 const rvec* p1 = s1->cg_p.rvec_array();
709 rvec* p2 = s2->cg_p.rvec_array();
710 #pragma omp for schedule(static) nowait
711 for (int i = start; i < end; i++)
713 // Trivial OpenMP block that does not throw
714 copy_rvec(p1[i], p2[i]);
718 if (DOMAINDECOMP(cr))
720 /* OpenMP does not supported unsigned loop variables */
721 #pragma omp for schedule(static) nowait
722 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
724 s2->cg_gl[i] = s1->cg_gl[i];
729 if (DOMAINDECOMP(cr))
731 s2->ddp_count = s1->ddp_count;
732 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
738 validStep = constr->apply(TRUE,
743 s1->x.arrayRefWithPadding(),
744 s2->x.arrayRefWithPadding(),
747 s2->lambda[FreeEnergyPerturbationCouplingType::Bonded],
749 gmx::ArrayRefWithPadding<RVec>(),
752 gmx::ConstraintVariable::Positions);
756 /* This global reduction will affect performance at high
757 * parallelization, but we can not really avoid it.
758 * But usually EM is not run at high parallelization.
760 int reductionBuffer = static_cast<int>(!validStep);
761 gmx_sumi(1, &reductionBuffer, cr);
762 validStep = (reductionBuffer == 0);
765 // We should move this check to the different minimizers
766 if (!validStep && ir->eI != IntegrationAlgorithm::Steep)
769 "The coordinates could not be constrained. Minimizer '%s' can not handle "
770 "constraint failures, use minimizer '%s' before using '%s'.",
771 enumValueToString(ir->eI),
772 enumValueToString(IntegrationAlgorithm::Steep),
773 enumValueToString(ir->eI));
780 //! Prepare EM for using domain decomposition parallellization
781 static void em_dd_partition_system(FILE* fplog,
782 const gmx::MDLogger& mdlog,
785 const gmx_mtop_t& top_global,
786 const t_inputrec* ir,
787 gmx::ImdSession* imdSession,
791 gmx::MDAtoms* mdAtoms,
793 VirtualSitesHandler* vsite,
794 gmx::Constraints* constr,
796 gmx_wallcycle_t wcycle)
798 /* Repartition the domain decomposition */
799 dd_partition_system(fplog,
820 dd_store_state(*cr->dd, &ems->s);
826 /*! \brief Class to handle the work of setting and doing an energy evaluation.
828 * This class is a mere aggregate of parameters to pass to evaluate an
829 * energy, so that future changes to names and types of them consume
830 * less time when refactoring other code.
832 * Aggregate initialization is used, for which the chief risk is that
833 * if a member is added at the end and not all initializer lists are
834 * updated, then the member will be value initialized, which will
835 * typically mean initialization to zero.
837 * Use a braced initializer list to construct one of these. */
838 class EnergyEvaluator
841 /*! \brief Evaluates an energy on the state in \c ems.
843 * \todo In practice, the same objects mu_tot, vir, and pres
844 * are always passed to this function, so we would rather have
845 * them as data members. However, their C-array types are
846 * unsuited for aggregate initialization. When the types
847 * improve, the call signature of this method can be reduced.
849 void run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst);
850 //! Handles logging (deprecated).
853 const gmx::MDLogger& mdlog;
854 //! Handles communication.
856 //! Coordinates multi-simulations.
857 const gmx_multisim_t* ms;
858 //! Holds the simulation topology.
859 const gmx_mtop_t& top_global;
860 //! Holds the domain topology.
862 //! User input options.
863 const t_inputrec* inputrec;
864 //! The Interactive Molecular Dynamics session.
865 gmx::ImdSession* imdSession;
866 //! The pull work object.
868 //! Manages flop accounting.
870 //! Manages wall cycle accounting.
871 gmx_wallcycle_t wcycle;
872 //! Coordinates global reduction.
873 gmx_global_stat_t gstat;
874 //! Handles virtual sites.
875 VirtualSitesHandler* vsite;
876 //! Handles constraints.
877 gmx::Constraints* constr;
878 //! Per-atom data for this domain.
879 gmx::MDAtoms* mdAtoms;
880 //! Handles how to calculate the forces.
882 //! Schedule of force-calculation work each step for this task.
883 MdrunScheduleWorkload* runScheduleWork;
884 //! Stores the computed energies.
885 gmx_enerdata_t* enerd;
888 void EnergyEvaluator::run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst)
892 tensor force_vir, shake_vir, ekin;
896 /* Set the time to the initial time, the time does not change during EM */
897 t = inputrec->init_t;
899 if (bFirst || (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
901 /* This is the first state or an old state used before the last ns */
907 if (inputrec->nstlist > 0)
915 vsite->construct(ems->s.x, {}, ems->s.box, gmx::VSiteOperation::Positions);
918 if (DOMAINDECOMP(cr) && bNS)
920 /* Repartition the domain decomposition */
921 em_dd_partition_system(
922 fplog, mdlog, count, cr, top_global, inputrec, imdSession, pull_work, ems, top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
925 /* Calc force & energy on new trial position */
926 /* do_force always puts the charge groups in the box and shifts again
927 * We do not unshift, so molecules are always whole in congrad.c
942 ems->s.x.arrayRefWithPadding(),
955 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY
956 | (bNS ? GMX_FORCE_NS : 0),
957 DDBalanceRegionHandler(cr));
959 /* Clear the unused shake virial and pressure */
960 clear_mat(shake_vir);
963 /* Communicate stuff when parallel */
964 if (PAR(cr) && inputrec->eI != IntegrationAlgorithm::NM)
966 wallcycle_start(wcycle, ewcMoveE);
975 gmx::ArrayRef<real>{},
977 std::vector<real>(1, terminate),
979 CGLO_ENERGY | CGLO_PRESSURE | CGLO_CONSTRAINT);
981 wallcycle_stop(wcycle, ewcMoveE);
984 if (fr->dispersionCorrection)
986 /* Calculate long range corrections to pressure and energy */
987 const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
988 ems->s.box, ems->s.lambda[FreeEnergyPerturbationCouplingType::Vdw]);
990 enerd->term[F_DISPCORR] = correction.energy;
991 enerd->term[F_EPOT] += correction.energy;
992 enerd->term[F_PRES] += correction.pressure;
993 enerd->term[F_DVDL] += correction.dvdl;
997 enerd->term[F_DISPCORR] = 0;
1000 ems->epot = enerd->term[F_EPOT];
1004 /* Project out the constraint components of the force */
1005 bool needsLogging = false;
1006 bool computeEnergy = false;
1007 bool computeVirial = true;
1009 auto f = ems->f.view().forceWithPadding();
1010 constr->apply(needsLogging,
1015 ems->s.x.arrayRefWithPadding(),
1017 f.unpaddedArrayRef(),
1019 ems->s.lambda[FreeEnergyPerturbationCouplingType::Bonded],
1021 gmx::ArrayRefWithPadding<RVec>(),
1024 gmx::ConstraintVariable::ForceDispl);
1025 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1026 m_add(force_vir, shake_vir, vir);
1030 copy_mat(force_vir, vir);
1034 enerd->term[F_PRES] = calc_pres(fr->pbcType, inputrec->nwall, ems->s.box, ekin, vir, pres);
1036 if (inputrec->efep != FreeEnergyPerturbationType::No)
1038 accumulateKineticLambdaComponents(enerd, ems->s.lambda, *inputrec->fepvals);
1041 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
1043 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
1049 //! Parallel utility summing energies and forces
1050 static double reorder_partsum(const t_commrec* cr,
1051 const t_grpopts* opts,
1052 const gmx_mtop_t& top_global,
1053 const em_state_t* s_min,
1054 const em_state_t* s_b)
1058 fprintf(debug, "Doing reorder_partsum\n");
1061 auto fm = s_min->f.view().force();
1062 auto fb = s_b->f.view().force();
1064 /* Collect fm in a global vector fmg.
1065 * This conflicts with the spirit of domain decomposition,
1066 * but to fully optimize this a much more complicated algorithm is required.
1068 const int natoms = top_global.natoms;
1072 gmx::ArrayRef<const int> indicesMin = s_min->s.cg_gl;
1074 for (int a : indicesMin)
1076 copy_rvec(fm[i], fmg[a]);
1079 gmx_sum(top_global.natoms * 3, fmg[0], cr);
1081 /* Now we will determine the part of the sum for the cgs in state s_b */
1082 gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
1087 gmx::ArrayRef<const unsigned char> grpnrFREEZE =
1088 top_global.groups.groupNumbers[SimulationAtomGroupType::Freeze];
1089 for (int a : indicesB)
1091 if (!grpnrFREEZE.empty())
1093 gf = grpnrFREEZE[i];
1095 for (int m = 0; m < DIM; m++)
1097 if (!opts->nFreeze[gf][m])
1099 partsum += (fb[i][m] - fmg[a][m]) * fb[i][m];
1110 //! Print some stuff, like beta, whatever that means.
1111 static real pr_beta(const t_commrec* cr,
1112 const t_grpopts* opts,
1114 const gmx_mtop_t& top_global,
1115 const em_state_t* s_min,
1116 const em_state_t* s_b)
1120 /* This is just the classical Polak-Ribiere calculation of beta;
1121 * it looks a bit complicated since we take freeze groups into account,
1122 * and might have to sum it in parallel runs.
1125 if (!DOMAINDECOMP(cr)
1126 || (s_min->s.ddp_count == cr->dd->ddp_count && s_b->s.ddp_count == cr->dd->ddp_count))
1128 auto fm = s_min->f.view().force();
1129 auto fb = s_b->f.view().force();
1132 /* This part of code can be incorrect with DD,
1133 * since the atom ordering in s_b and s_min might differ.
1135 for (int i = 0; i < mdatoms->homenr; i++)
1137 if (mdatoms->cFREEZE)
1139 gf = mdatoms->cFREEZE[i];
1141 for (int m = 0; m < DIM; m++)
1143 if (!opts->nFreeze[gf][m])
1145 sum += (fb[i][m] - fm[i][m]) * fb[i][m];
1152 /* We need to reorder cgs while summing */
1153 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1157 gmx_sumd(1, &sum, cr);
1160 return sum / gmx::square(s_min->fnorm);
1166 void LegacySimulator::do_cg()
1168 const char* CG = "Polak-Ribiere Conjugate Gradients";
1170 gmx_localtop_t top(top_global.ffparams);
1171 gmx_global_stat_t gstat;
1172 double tmp, minstep;
1174 real a, b, c, beta = 0.0;
1177 gmx_bool converged, foundlower;
1178 rvec mu_tot = { 0 };
1179 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1181 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1182 int m, step, nminstep;
1183 auto mdatoms = mdAtoms->mdatoms();
1188 "Note that activating conjugate gradient energy minimization via the "
1189 "integrator .mdp option and the command gmx mdrun may "
1190 "be available in a different form in a future version of GROMACS, "
1191 "e.g. gmx minimize and an .mdp option.");
1197 // In CG, the state is extended with a search direction
1198 state_global->flags |= enumValueToBitMask(StateEntry::Cgp);
1200 // Ensure the extra per-atom state array gets allocated
1201 state_change_natoms(state_global, state_global->natoms);
1203 // Initialize the search direction to zero
1204 for (RVec& cg_p : state_global->cg_p)
1210 /* Create 4 states on the stack and extract pointers that we will swap */
1211 em_state_t s0{}, s1{}, s2{}, s3{};
1212 em_state_t* s_min = &s0;
1213 em_state_t* s_a = &s1;
1214 em_state_t* s_b = &s2;
1215 em_state_t* s_c = &s3;
1217 /* Init em and store the local state in s_min */
1236 const bool simulationsShareState = false;
1237 gmx_mdoutf* outf = init_mdoutf(fplog,
1248 StartingBehavior::NewSimulation,
1249 simulationsShareState,
1251 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
1257 StartingBehavior::NewSimulation,
1258 simulationsShareState,
1261 /* Print to log file */
1262 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1264 /* Max number of steps */
1265 number_steps = inputrec->nsteps;
1269 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1273 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1276 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
1277 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1278 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
1279 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1280 /* do_force always puts the charge groups in the box and shifts again
1281 * We do not unshift, so molecules are always whole in congrad.c
1283 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1287 /* Copy stuff to the energy bin for easy printing etc. */
1288 matrix nullBox = {};
1289 energyOutput.addDataAtEnergyStep(false,
1291 static_cast<double>(step),
1307 EnergyOutput::printHeader(fplog, step, step);
1308 energyOutput.printStepToEnergyFile(
1309 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
1312 /* Estimate/guess the initial stepsize */
1313 stepsize = inputrec->em_stepsize / s_min->fnorm;
1317 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1318 fprintf(stderr, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1319 fprintf(stderr, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1320 fprintf(stderr, "\n");
1321 /* and copy to the log file too... */
1322 fprintf(fplog, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1323 fprintf(fplog, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1324 fprintf(fplog, "\n");
1326 /* Start the loop over CG steps.
1327 * Each successful step is counted, and we continue until
1328 * we either converge or reach the max number of steps.
1331 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1334 /* start taking steps in a new direction
1335 * First time we enter the routine, beta=0, and the direction is
1336 * simply the negative gradient.
1339 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1340 gmx::ArrayRef<gmx::RVec> pm = s_min->s.cg_p;
1341 gmx::ArrayRef<const gmx::RVec> sfm = s_min->f.view().force();
1344 for (int i = 0; i < mdatoms->homenr; i++)
1346 if (mdatoms->cFREEZE)
1348 gf = mdatoms->cFREEZE[i];
1350 for (m = 0; m < DIM; m++)
1352 if (!inputrec->opts.nFreeze[gf][m])
1354 pm[i][m] = sfm[i][m] + beta * pm[i][m];
1355 gpa -= pm[i][m] * sfm[i][m];
1356 /* f is negative gradient, thus the sign */
1365 /* Sum the gradient along the line across CPUs */
1368 gmx_sumd(1, &gpa, cr);
1371 /* Calculate the norm of the search vector */
1372 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1374 /* Just in case stepsize reaches zero due to numerical precision... */
1377 stepsize = inputrec->em_stepsize / pnorm;
1381 * Double check the value of the derivative in the search direction.
1382 * If it is positive it must be due to the old information in the
1383 * CG formula, so just remove that and start over with beta=0.
1384 * This corresponds to a steepest descent step.
1389 step--; /* Don't count this step since we are restarting */
1390 continue; /* Go back to the beginning of the big for-loop */
1393 /* Calculate minimum allowed stepsize, before the average (norm)
1394 * relative change in coordinate is smaller than precision
1397 auto s_min_x = makeArrayRef(s_min->s.x);
1398 for (int i = 0; i < mdatoms->homenr; i++)
1400 for (m = 0; m < DIM; m++)
1402 tmp = fabs(s_min_x[i][m]);
1407 tmp = pm[i][m] / tmp;
1408 minstep += tmp * tmp;
1411 /* Add up from all CPUs */
1414 gmx_sumd(1, &minstep, cr);
1417 minstep = GMX_REAL_EPS / sqrt(minstep / (3 * top_global.natoms));
1419 if (stepsize < minstep)
1425 /* Write coordinates if necessary */
1426 do_x = do_per_step(step, inputrec->nstxout);
1427 do_f = do_per_step(step, inputrec->nstfout);
1430 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, step, s_min, state_global, observablesHistory);
1432 /* Take a step downhill.
1433 * In theory, we should minimize the function along this direction.
1434 * That is quite possible, but it turns out to take 5-10 function evaluations
1435 * for each line. However, we dont really need to find the exact minimum -
1436 * it is much better to start a new CG step in a modified direction as soon
1437 * as we are close to it. This will save a lot of energy evaluations.
1439 * In practice, we just try to take a single step.
1440 * If it worked (i.e. lowered the energy), we increase the stepsize but
1441 * the continue straight to the next CG step without trying to find any minimum.
1442 * If it didn't work (higher energy), there must be a minimum somewhere between
1443 * the old position and the new one.
1445 * Due to the finite numerical accuracy, it turns out that it is a good idea
1446 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1447 * This leads to lower final energies in the tests I've done. / Erik
1449 s_a->epot = s_min->epot;
1451 c = a + stepsize; /* reference position along line is zero */
1453 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1455 em_dd_partition_system(fplog,
1473 /* Take a trial step (new coords in s_c) */
1474 do_em_step(cr, inputrec, mdatoms, s_min, c, s_min->s.cg_p.constArrayRefWithPadding(), s_c, constr, -1);
1477 /* Calculate energy for the trial step */
1478 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1480 /* Calc derivative along line */
1481 const rvec* pc = s_c->s.cg_p.rvec_array();
1482 gmx::ArrayRef<const gmx::RVec> sfc = s_c->f.view().force();
1484 for (int i = 0; i < mdatoms->homenr; i++)
1486 for (m = 0; m < DIM; m++)
1488 gpc -= pc[i][m] * sfc[i][m]; /* f is negative gradient, thus the sign */
1491 /* Sum the gradient along the line across CPUs */
1494 gmx_sumd(1, &gpc, cr);
1497 /* This is the max amount of increase in energy we tolerate */
1498 tmp = std::sqrt(GMX_REAL_EPS) * fabs(s_a->epot);
1500 /* Accept the step if the energy is lower, or if it is not significantly higher
1501 * and the line derivative is still negative.
1503 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1506 /* Great, we found a better energy. Increase step for next iteration
1507 * if we are still going down, decrease it otherwise
1511 stepsize *= 1.618034; /* The golden section */
1515 stepsize *= 0.618034; /* 1/golden section */
1520 /* New energy is the same or higher. We will have to do some work
1521 * to find a smaller value in the interval. Take smaller step next time!
1524 stepsize *= 0.618034;
1528 /* OK, if we didn't find a lower value we will have to locate one now - there must
1529 * be one in the interval [a=0,c].
1530 * The same thing is valid here, though: Don't spend dozens of iterations to find
1531 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1532 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1534 * I also have a safeguard for potentially really pathological functions so we never
1535 * take more than 20 steps before we give up ...
1537 * If we already found a lower value we just skip this step and continue to the update.
1546 /* Select a new trial point.
1547 * If the derivatives at points a & c have different sign we interpolate to zero,
1548 * otherwise just do a bisection.
1550 if (gpa < 0 && gpc > 0)
1552 b = a + gpa * (a - c) / (gpc - gpa);
1559 /* safeguard if interpolation close to machine accuracy causes errors:
1560 * never go outside the interval
1562 if (b <= a || b >= c)
1567 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1569 /* Reload the old state */
1570 em_dd_partition_system(fplog,
1588 /* Take a trial step to this new point - new coords in s_b */
1589 do_em_step(cr, inputrec, mdatoms, s_min, b, s_min->s.cg_p.constArrayRefWithPadding(), s_b, constr, -1);
1592 /* Calculate energy for the trial step */
1593 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1595 /* p does not change within a step, but since the domain decomposition
1596 * might change, we have to use cg_p of s_b here.
1598 const rvec* pb = s_b->s.cg_p.rvec_array();
1599 gmx::ArrayRef<const gmx::RVec> sfb = s_b->f.view().force();
1601 for (int i = 0; i < mdatoms->homenr; i++)
1603 for (m = 0; m < DIM; m++)
1605 gpb -= pb[i][m] * sfb[i][m]; /* f is negative gradient, thus the sign */
1608 /* Sum the gradient along the line across CPUs */
1611 gmx_sumd(1, &gpb, cr);
1616 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a->epot, s_b->epot, s_c->epot, gpb);
1619 epot_repl = s_b->epot;
1621 /* Keep one of the intervals based on the value of the derivative at the new point */
1624 /* Replace c endpoint with b */
1625 swap_em_state(&s_b, &s_c);
1631 /* Replace a endpoint with b */
1632 swap_em_state(&s_b, &s_a);
1638 * Stop search as soon as we find a value smaller than the endpoints.
1639 * Never run more than 20 steps, no matter what.
1642 } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) && (nminstep < 20));
1644 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot) * GMX_REAL_EPS || nminstep >= 20)
1646 /* OK. We couldn't find a significantly lower energy.
1647 * If beta==0 this was steepest descent, and then we give up.
1648 * If not, set beta=0 and restart with steepest descent before quitting.
1658 /* Reset memory before giving up */
1664 /* Select min energy state of A & C, put the best in B.
1666 if (s_c->epot < s_a->epot)
1670 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c->epot, s_a->epot);
1672 swap_em_state(&s_b, &s_c);
1679 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a->epot, s_c->epot);
1681 swap_em_state(&s_b, &s_a);
1689 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n", s_c->epot);
1691 swap_em_state(&s_b, &s_c);
1695 /* new search direction */
1696 /* beta = 0 means forget all memory and restart with steepest descents. */
1697 if (nstcg && ((step % nstcg) == 0))
1703 /* s_min->fnorm cannot be zero, because then we would have converged
1707 /* Polak-Ribiere update.
1708 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1710 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1712 /* Limit beta to prevent oscillations */
1713 if (fabs(beta) > 5.0)
1719 /* update positions */
1720 swap_em_state(&s_min, &s_b);
1723 /* Print it if necessary */
1726 if (mdrunOptions.verbose)
1728 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1730 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1733 s_min->fnorm / sqrtNumAtoms,
1738 /* Store the new (lower) energies */
1739 matrix nullBox = {};
1740 energyOutput.addDataAtEnergyStep(false,
1742 static_cast<double>(step),
1758 do_log = do_per_step(step, inputrec->nstlog);
1759 do_ene = do_per_step(step, inputrec->nstenergy);
1761 imdSession->fillEnergyRecord(step, TRUE);
1765 EnergyOutput::printHeader(fplog, step, step);
1767 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1771 do_log ? fplog : nullptr,
1778 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1779 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x, 0))
1781 imdSession->sendPositionsAndEnergies();
1784 /* Stop when the maximum force lies below tolerance.
1785 * If we have reached machine precision, converged is already set to true.
1787 converged = converged || (s_min->fmax < inputrec->em_tol);
1789 } /* End of the loop */
1793 step--; /* we never took that last step in this case */
1795 if (s_min->fmax > inputrec->em_tol)
1799 warn_step(fplog, inputrec->em_tol, s_min->fmax, step - 1 == number_steps, FALSE);
1806 /* If we printed energy and/or logfile last step (which was the last step)
1807 * we don't have to do it again, but otherwise print the final values.
1811 /* Write final value to log since we didn't do anything the last step */
1812 EnergyOutput::printHeader(fplog, step, step);
1814 if (!do_ene || !do_log)
1816 /* Write final energy file entries */
1817 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1821 !do_log ? fplog : nullptr,
1829 /* Print some stuff... */
1832 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1836 * For accurate normal mode calculation it is imperative that we
1837 * store the last conformation into the full precision binary trajectory.
1839 * However, we should only do it if we did NOT already write this step
1840 * above (which we did if do_x or do_f was true).
1842 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1843 * in the trajectory with the same step number.
1845 do_x = !do_per_step(step, inputrec->nstxout);
1846 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1849 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, s_min, state_global, observablesHistory);
1854 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1855 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1856 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1858 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1861 finish_em(cr, outf, walltime_accounting, wcycle);
1863 /* To print the actual number of steps we needed somewhere */
1864 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1868 void LegacySimulator::do_lbfgs()
1870 static const char* LBFGS = "Low-Memory BFGS Minimizer";
1872 gmx_localtop_t top(top_global.ffparams);
1873 gmx_global_stat_t gstat;
1874 auto mdatoms = mdAtoms->mdatoms();
1879 "Note that activating L-BFGS energy minimization via the "
1880 "integrator .mdp option and the command gmx mdrun may "
1881 "be available in a different form in a future version of GROMACS, "
1882 "e.g. gmx minimize and an .mdp option.");
1886 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1889 if (nullptr != constr)
1893 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1894 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1897 const int n = 3 * state_global->natoms;
1898 const int nmaxcorr = inputrec->nbfgscorr;
1900 std::vector<real> p(n);
1901 std::vector<real> rho(nmaxcorr);
1902 std::vector<real> alpha(nmaxcorr);
1904 std::vector<std::vector<real>> dx(nmaxcorr);
1905 for (auto& dxCorr : dx)
1910 std::vector<std::vector<real>> dg(nmaxcorr);
1911 for (auto& dgCorr : dg)
1938 const bool simulationsShareState = false;
1939 gmx_mdoutf* outf = init_mdoutf(fplog,
1950 StartingBehavior::NewSimulation,
1951 simulationsShareState,
1953 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
1959 StartingBehavior::NewSimulation,
1960 simulationsShareState,
1963 const int start = 0;
1964 const int end = mdatoms->homenr;
1966 /* We need 4 working states */
1967 em_state_t s0{}, s1{}, s2{}, s3{};
1968 em_state_t* sa = &s0;
1969 em_state_t* sb = &s1;
1970 em_state_t* sc = &s2;
1971 em_state_t* last = &s3;
1972 /* Initialize by copying the state from ems (we could skip x and f here) */
1977 /* Print to log file */
1978 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1980 /* Max number of steps */
1981 const int number_steps = inputrec->nsteps;
1983 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1984 std::vector<bool> frozen(n);
1986 for (int i = start; i < end; i++)
1988 if (mdatoms->cFREEZE)
1990 gf = mdatoms->cFREEZE[i];
1992 for (int m = 0; m < DIM; m++)
1994 frozen[3 * i + m] = (inputrec->opts.nFreeze[gf][m] != 0);
1999 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
2003 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
2008 vsite->construct(state_global->x, {}, state_global->box, VSiteOperation::Positions);
2011 /* Call the force routine and some auxiliary (neighboursearching etc.) */
2012 /* do_force always puts the charge groups in the box and shifts again
2013 * We do not unshift, so molecules are always whole
2016 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2017 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2018 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2022 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
2026 /* Copy stuff to the energy bin for easy printing etc. */
2027 matrix nullBox = {};
2028 energyOutput.addDataAtEnergyStep(false,
2030 static_cast<double>(step),
2046 EnergyOutput::printHeader(fplog, step, step);
2047 energyOutput.printStepToEnergyFile(
2048 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
2051 /* Set the initial step.
2052 * since it will be multiplied by the non-normalized search direction
2053 * vector (force vector the first time), we scale it by the
2054 * norm of the force.
2059 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2060 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2061 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2062 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2063 fprintf(stderr, "\n");
2064 /* and copy to the log file too... */
2065 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2066 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2067 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2068 fprintf(fplog, "\n");
2071 // Point is an index to the memory of search directions, where 0 is the first one.
2074 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
2075 real* fInit = static_cast<real*>(ems.f.view().force().data()[0]);
2076 for (int i = 0; i < n; i++)
2080 dx[point][i] = fInit[i]; /* Initial search direction */
2088 // Stepsize will be modified during the search, and actually it is not critical
2089 // (the main efficiency in the algorithm comes from changing directions), but
2090 // we still need an initial value, so estimate it as the inverse of the norm
2091 // so we take small steps where the potential fluctuates a lot.
2092 double stepsize = 1.0 / ems.fnorm;
2094 /* Start the loop over BFGS steps.
2095 * Each successful step is counted, and we continue until
2096 * we either converge or reach the max number of steps.
2104 /* Set the gradient from the force */
2105 bool converged = false;
2106 for (int step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
2109 /* Write coordinates if necessary */
2110 const bool do_x = do_per_step(step, inputrec->nstxout);
2111 const bool do_f = do_per_step(step, inputrec->nstfout);
2116 mdof_flags |= MDOF_X;
2121 mdof_flags |= MDOF_F;
2126 mdof_flags |= MDOF_IMD;
2129 gmx::WriteCheckpointDataHolder checkpointDataHolder;
2130 mdoutf_write_to_trajectory_files(fplog,
2136 static_cast<real>(step),
2140 ems.f.view().force(),
2141 &checkpointDataHolder);
2143 /* Do the linesearching in the direction dx[point][0..(n-1)] */
2145 /* make s a pointer to current search direction - point=0 first time we get here */
2146 gmx::ArrayRef<const real> s = dx[point];
2148 const real* xx = static_cast<real*>(ems.s.x.rvec_array()[0]);
2149 const real* ff = static_cast<real*>(ems.f.view().force().data()[0]);
2151 // calculate line gradient in position A
2153 for (int i = 0; i < n; i++)
2155 gpa -= s[i] * ff[i];
2158 /* Calculate minimum allowed stepsize along the line, before the average (norm)
2159 * relative change in coordinate is smaller than precision
2162 for (int i = 0; i < n; i++)
2164 double tmp = fabs(xx[i]);
2170 minstep += tmp * tmp;
2172 minstep = GMX_REAL_EPS / sqrt(minstep / n);
2174 if (stepsize < minstep)
2180 // Before taking any steps along the line, store the old position
2182 real* lastx = static_cast<real*>(last->s.x.data()[0]);
2183 real* lastf = static_cast<real*>(last->f.view().force().data()[0]);
2184 const real Epot0 = ems.epot;
2188 /* Take a step downhill.
2189 * In theory, we should find the actual minimum of the function in this
2190 * direction, somewhere along the line.
2191 * That is quite possible, but it turns out to take 5-10 function evaluations
2192 * for each line. However, we dont really need to find the exact minimum -
2193 * it is much better to start a new BFGS step in a modified direction as soon
2194 * as we are close to it. This will save a lot of energy evaluations.
2196 * In practice, we just try to take a single step.
2197 * If it worked (i.e. lowered the energy), we increase the stepsize but
2198 * continue straight to the next BFGS step without trying to find any minimum,
2199 * i.e. we change the search direction too. If the line was smooth, it is
2200 * likely we are in a smooth region, and then it makes sense to take longer
2201 * steps in the modified search direction too.
2203 * If it didn't work (higher energy), there must be a minimum somewhere between
2204 * the old position and the new one. Then we need to start by finding a lower
2205 * value before we change search direction. Since the energy was apparently
2206 * quite rough, we need to decrease the step size.
2208 * Due to the finite numerical accuracy, it turns out that it is a good idea
2209 * to accept a SMALL increase in energy, if the derivative is still downhill.
2210 * This leads to lower final energies in the tests I've done. / Erik
2213 // State "A" is the first position along the line.
2214 // reference position along line is initially zero
2217 // Check stepsize first. We do not allow displacements
2218 // larger than emstep.
2224 // Pick a new position C by adding stepsize to A.
2227 // Calculate what the largest change in any individual coordinate
2228 // would be (translation along line * gradient along line)
2230 for (int i = 0; i < n; i++)
2232 real delta = c * s[i];
2233 if (delta > maxdelta)
2238 // If any displacement is larger than the stepsize limit, reduce the step
2239 if (maxdelta > inputrec->em_stepsize)
2243 } while (maxdelta > inputrec->em_stepsize);
2245 // Take a trial step and move the coordinate array xc[] to position C
2246 real* xc = static_cast<real*>(sc->s.x.rvec_array()[0]);
2247 for (int i = 0; i < n; i++)
2249 xc[i] = lastx[i] + c * s[i];
2253 // Calculate energy for the trial step in position C
2254 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2256 // Calc line gradient in position C
2257 real* fc = static_cast<real*>(sc->f.view().force()[0]);
2259 for (int i = 0; i < n; i++)
2261 gpc -= s[i] * fc[i]; /* f is negative gradient, thus the sign */
2263 /* Sum the gradient along the line across CPUs */
2266 gmx_sumd(1, &gpc, cr);
2269 // This is the max amount of increase in energy we tolerate.
2270 // By allowing VERY small changes (close to numerical precision) we
2271 // frequently find even better (lower) final energies.
2272 double tmp = std::sqrt(GMX_REAL_EPS) * fabs(sa->epot);
2274 // Accept the step if the energy is lower in the new position C (compared to A),
2275 // or if it is not significantly higher and the line derivative is still negative.
2276 bool foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2277 // If true, great, we found a better energy. We no longer try to alter the
2278 // stepsize, but simply accept this new better position. The we select a new
2279 // search direction instead, which will be much more efficient than continuing
2280 // to take smaller steps along a line. Set fnorm based on the new C position,
2281 // which will be used to update the stepsize to 1/fnorm further down.
2283 // If false, the energy is NOT lower in point C, i.e. it will be the same
2284 // or higher than in point A. In this case it is pointless to move to point C,
2285 // so we will have to do more iterations along the same line to find a smaller
2286 // value in the interval [A=0.0,C].
2287 // Here, A is still 0.0, but that will change when we do a search in the interval
2288 // [0.0,C] below. That search we will do by interpolation or bisection rather
2289 // than with the stepsize, so no need to modify it. For the next search direction
2290 // it will be reset to 1/fnorm anyway.
2295 // OK, if we didn't find a lower value we will have to locate one now - there must
2296 // be one in the interval [a,c].
2297 // The same thing is valid here, though: Don't spend dozens of iterations to find
2298 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2299 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2300 // I also have a safeguard for potentially really pathological functions so we never
2301 // take more than 20 steps before we give up.
2302 // If we already found a lower value we just skip this step and continue to the update.
2307 // Select a new trial point B in the interval [A,C].
2308 // If the derivatives at points a & c have different sign we interpolate to zero,
2309 // otherwise just do a bisection since there might be multiple minima/maxima
2310 // inside the interval.
2312 if (gpa < 0 && gpc > 0)
2314 b = a + gpa * (a - c) / (gpc - gpa);
2321 /* safeguard if interpolation close to machine accuracy causes errors:
2322 * never go outside the interval
2324 if (b <= a || b >= c)
2329 // Take a trial step to point B
2330 real* xb = static_cast<real*>(sb->s.x.rvec_array()[0]);
2331 for (int i = 0; i < n; i++)
2333 xb[i] = lastx[i] + b * s[i];
2337 // Calculate energy for the trial step in point B
2338 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2341 // Calculate gradient in point B
2342 real* fb = static_cast<real*>(sb->f.view().force()[0]);
2344 for (int i = 0; i < n; i++)
2346 gpb -= s[i] * fb[i]; /* f is negative gradient, thus the sign */
2348 /* Sum the gradient along the line across CPUs */
2351 gmx_sumd(1, &gpb, cr);
2354 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2355 // at the new point B, and rename the endpoints of this new interval A and C.
2358 /* Replace c endpoint with b */
2360 /* copy state b to c */
2365 /* Replace a endpoint with b */
2367 /* copy state b to a */
2372 * Stop search as soon as we find a value smaller than the endpoints,
2373 * or if the tolerance is below machine precision.
2374 * Never run more than 20 steps, no matter what.
2377 } while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2379 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2381 /* OK. We couldn't find a significantly lower energy.
2382 * If ncorr==0 this was steepest descent, and then we give up.
2383 * If not, reset memory to restart as steepest descent before quitting.
2395 /* Search in gradient direction */
2396 for (int i = 0; i < n; i++)
2398 dx[point][i] = ff[i];
2400 /* Reset stepsize */
2401 stepsize = 1.0 / fnorm;
2406 /* Select min energy state of A & C, put the best in xx/ff/Epot
2408 if (sc->epot < sa->epot)
2429 /* Update the memory information, and calculate a new
2430 * approximation of the inverse hessian
2433 /* Have new data in Epot, xx, ff */
2434 if (ncorr < nmaxcorr)
2439 for (int i = 0; i < n; i++)
2441 dg[point][i] = lastf[i] - ff[i];
2442 dx[point][i] *= step_taken;
2447 for (int i = 0; i < n; i++)
2449 dgdg += dg[point][i] * dg[point][i];
2450 dgdx += dg[point][i] * dx[point][i];
2453 const real diag = dgdx / dgdg;
2455 rho[point] = 1.0 / dgdx;
2458 if (point >= nmaxcorr)
2464 for (int i = 0; i < n; i++)
2471 /* Recursive update. First go back over the memory points */
2472 for (int k = 0; k < ncorr; k++)
2481 for (int i = 0; i < n; i++)
2483 sq += dx[cp][i] * p[i];
2486 alpha[cp] = rho[cp] * sq;
2488 for (int i = 0; i < n; i++)
2490 p[i] -= alpha[cp] * dg[cp][i];
2494 for (int i = 0; i < n; i++)
2499 /* And then go forward again */
2500 for (int k = 0; k < ncorr; k++)
2503 for (int i = 0; i < n; i++)
2505 yr += p[i] * dg[cp][i];
2508 real beta = rho[cp] * yr;
2509 beta = alpha[cp] - beta;
2511 for (int i = 0; i < n; i++)
2513 p[i] += beta * dx[cp][i];
2523 for (int i = 0; i < n; i++)
2527 dx[point][i] = p[i];
2535 /* Print it if necessary */
2538 if (mdrunOptions.verbose)
2540 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2542 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2545 ems.fnorm / sqrtNumAtoms,
2550 /* Store the new (lower) energies */
2551 matrix nullBox = {};
2552 energyOutput.addDataAtEnergyStep(false,
2554 static_cast<double>(step),
2570 do_log = do_per_step(step, inputrec->nstlog);
2571 do_ene = do_per_step(step, inputrec->nstenergy);
2573 imdSession->fillEnergyRecord(step, TRUE);
2577 EnergyOutput::printHeader(fplog, step, step);
2579 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2583 do_log ? fplog : nullptr,
2590 /* Send x and E to IMD client, if bIMD is TRUE. */
2591 if (imdSession->run(step, TRUE, state_global->box, state_global->x, 0) && MASTER(cr))
2593 imdSession->sendPositionsAndEnergies();
2596 // Reset stepsize in we are doing more iterations
2599 /* Stop when the maximum force lies below tolerance.
2600 * If we have reached machine precision, converged is already set to true.
2602 converged = converged || (ems.fmax < inputrec->em_tol);
2604 } /* End of the loop */
2608 step--; /* we never took that last step in this case */
2610 if (ems.fmax > inputrec->em_tol)
2614 warn_step(fplog, inputrec->em_tol, ems.fmax, step - 1 == number_steps, FALSE);
2619 /* If we printed energy and/or logfile last step (which was the last step)
2620 * we don't have to do it again, but otherwise print the final values.
2622 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2624 EnergyOutput::printHeader(fplog, step, step);
2626 if (!do_ene || !do_log) /* Write final energy file entries */
2628 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2632 !do_log ? fplog : nullptr,
2639 /* Print some stuff... */
2642 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2646 * For accurate normal mode calculation it is imperative that we
2647 * store the last conformation into the full precision binary trajectory.
2649 * However, we should only do it if we did NOT already write this step
2650 * above (which we did if do_x or do_f was true).
2652 const bool do_x = !do_per_step(step, inputrec->nstxout);
2653 const bool do_f = !do_per_step(step, inputrec->nstfout);
2655 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, &ems, state_global, observablesHistory);
2659 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2660 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2661 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2663 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2666 finish_em(cr, outf, walltime_accounting, wcycle);
2668 /* To print the actual number of steps we needed somewhere */
2669 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2672 void LegacySimulator::do_steep()
2674 const char* SD = "Steepest Descents";
2675 gmx_localtop_t top(top_global.ffparams);
2676 gmx_global_stat_t gstat;
2679 gmx_bool bDone, bAbort, do_x, do_f;
2681 rvec mu_tot = { 0 };
2684 int steps_accepted = 0;
2685 auto mdatoms = mdAtoms->mdatoms();
2690 "Note that activating steepest-descent energy minimization via the "
2691 "integrator .mdp option and the command gmx mdrun may "
2692 "be available in a different form in a future version of GROMACS, "
2693 "e.g. gmx minimize and an .mdp option.");
2695 /* Create 2 states on the stack and extract pointers that we will swap */
2696 em_state_t s0{}, s1{};
2697 em_state_t* s_min = &s0;
2698 em_state_t* s_try = &s1;
2700 /* Init em and store the local state in s_try */
2719 const bool simulationsShareState = false;
2720 gmx_mdoutf* outf = init_mdoutf(fplog,
2731 StartingBehavior::NewSimulation,
2732 simulationsShareState,
2734 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
2740 StartingBehavior::NewSimulation,
2741 simulationsShareState,
2744 /* Print to log file */
2745 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2747 /* Set variables for stepsize (in nm). This is the largest
2748 * step that we are going to make in any direction.
2750 ustep = inputrec->em_stepsize;
2753 /* Max number of steps */
2754 nsteps = inputrec->nsteps;
2758 /* Print to the screen */
2759 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2763 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2765 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2766 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2767 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2769 /**** HERE STARTS THE LOOP ****
2770 * count is the counter for the number of steps
2771 * bDone will be TRUE when the minimization has converged
2772 * bAbort will be TRUE when nsteps steps have been performed or when
2773 * the stepsize becomes smaller than is reasonable for machine precision
2778 while (!bDone && !bAbort)
2780 bAbort = (nsteps >= 0) && (count == nsteps);
2782 /* set new coordinates, except for first step */
2783 bool validStep = true;
2786 validStep = do_em_step(
2787 cr, inputrec, mdatoms, s_min, stepsize, s_min->f.view().forceWithPadding(), s_try, constr, count);
2792 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2796 // Signal constraint error during stepping with energy=inf
2797 s_try->epot = std::numeric_limits<real>::infinity();
2802 EnergyOutput::printHeader(fplog, count, count);
2807 s_min->epot = s_try->epot;
2810 /* Print it if necessary */
2813 if (mdrunOptions.verbose)
2816 "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2822 ((count == 0) || (s_try->epot < s_min->epot)) ? '\n' : '\r');
2826 if ((count == 0) || (s_try->epot < s_min->epot))
2828 /* Store the new (lower) energies */
2829 matrix nullBox = {};
2830 energyOutput.addDataAtEnergyStep(false,
2832 static_cast<double>(count),
2848 imdSession->fillEnergyRecord(count, TRUE);
2850 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2851 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2852 energyOutput.printStepToEnergyFile(
2853 mdoutf_get_fp_ene(outf), TRUE, do_dr, do_or, fplog, count, count, fr->fcdata.get(), nullptr);
2858 /* Now if the new energy is smaller than the previous...
2859 * or if this is the first step!
2860 * or if we did random steps!
2863 if ((count == 0) || (s_try->epot < s_min->epot))
2867 /* Test whether the convergence criterion is met... */
2868 bDone = (s_try->fmax < inputrec->em_tol);
2870 /* Copy the arrays for force, positions and energy */
2871 /* The 'Min' array always holds the coords and forces of the minimal
2873 swap_em_state(&s_min, &s_try);
2879 /* Write to trn, if necessary */
2880 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2881 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2883 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, count, s_min, state_global, observablesHistory);
2887 /* If energy is not smaller make the step smaller... */
2890 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2892 /* Reload the old state */
2893 em_dd_partition_system(fplog,
2912 // If the force is very small after finishing minimization,
2913 // we risk dividing by zero when calculating the step size.
2914 // So we check first if the minimization has stopped before
2915 // trying to obtain a new step size.
2918 /* Determine new step */
2919 stepsize = ustep / s_min->fmax;
2922 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2924 if (count == nsteps || ustep < 1e-12)
2926 if (count == nsteps || ustep < 1e-6)
2931 warn_step(fplog, inputrec->em_tol, s_min->fmax, count == nsteps, constr != nullptr);
2936 /* Send IMD energies and positions, if bIMD is TRUE. */
2937 if (imdSession->run(count,
2939 MASTER(cr) ? state_global->box : nullptr,
2940 MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>(),
2944 imdSession->sendPositionsAndEnergies();
2948 } /* End of the loop */
2950 /* Print some data... */
2953 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2955 write_em_traj(fplog,
2959 inputrec->nstfout != 0,
2960 ftp2fn(efSTO, nfile, fnm),
2966 observablesHistory);
2970 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2972 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2973 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2976 finish_em(cr, outf, walltime_accounting, wcycle);
2978 /* To print the actual number of steps we needed somewhere */
2980 // TODO: Avoid changing inputrec (#3854)
2981 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
2982 nonConstInputrec->nsteps = count;
2985 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2988 void LegacySimulator::do_nm()
2990 const char* NM = "Normal Mode Analysis";
2992 gmx_localtop_t top(top_global.ffparams);
2993 gmx_global_stat_t gstat;
2995 rvec mu_tot = { 0 };
2997 gmx_bool bSparse; /* use sparse matrix storage format */
2999 gmx_sparsematrix_t* sparse_matrix = nullptr;
3000 real* full_matrix = nullptr;
3002 /* added with respect to mdrun */
3004 real der_range = 10.0 * std::sqrt(GMX_REAL_EPS);
3006 bool bIsMaster = MASTER(cr);
3007 auto mdatoms = mdAtoms->mdatoms();
3012 "Note that activating normal-mode analysis via the integrator "
3013 ".mdp option and the command gmx mdrun may "
3014 "be available in a different form in a future version of GROMACS, "
3015 "e.g. gmx normal-modes.");
3017 if (constr != nullptr)
3021 "Constraints present with Normal Mode Analysis, this combination is not supported");
3024 gmx_shellfc_t* shellfc;
3026 em_state_t state_work{};
3028 /* Init em and store the local state in state_minimum */
3047 const bool simulationsShareState = false;
3048 gmx_mdoutf* outf = init_mdoutf(fplog,
3059 StartingBehavior::NewSimulation,
3060 simulationsShareState,
3063 std::vector<int> atom_index = get_atom_index(top_global);
3064 std::vector<gmx::RVec> fneg(atom_index.size(), { 0, 0, 0 });
3065 snew(dfdx, atom_index.size());
3071 "NOTE: This version of GROMACS has been compiled in single precision,\n"
3072 " which MIGHT not be accurate enough for normal mode analysis.\n"
3073 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
3074 " are fairly modest even if you recompile in double precision.\n\n");
3078 /* Check if we can/should use sparse storage format.
3080 * Sparse format is only useful when the Hessian itself is sparse, which it
3081 * will be when we use a cutoff.
3082 * For small systems (n<1000) it is easier to always use full matrix format, though.
3084 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
3086 GMX_LOG(mdlog.warning)
3087 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
3090 else if (atom_index.size() < 1000)
3092 GMX_LOG(mdlog.warning)
3093 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
3099 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
3103 /* Number of dimensions, based on real atoms, that is not vsites or shell */
3104 sz = DIM * atom_index.size();
3106 fprintf(stderr, "Allocating Hessian memory...\n\n");
3110 sparse_matrix = gmx_sparsematrix_init(sz);
3111 sparse_matrix->compressed_symmetric = TRUE;
3115 snew(full_matrix, sz * sz);
3118 /* Write start time and temperature */
3119 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
3121 /* fudge nr of steps to nr of atoms */
3123 // TODO: Avoid changing inputrec (#3854)
3124 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
3125 nonConstInputrec->nsteps = atom_index.size() * 2;
3131 "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
3136 nnodes = cr->nnodes;
3138 /* Make evaluate_energy do a single node force calculation */
3140 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
3141 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
3142 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
3143 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
3144 cr->nnodes = nnodes;
3146 /* if forces are not small, warn user */
3147 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
3149 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
3150 if (state_work.fmax > 1.0e-3)
3152 GMX_LOG(mdlog.warning)
3154 "The force is probably not small enough to "
3155 "ensure that you are at a minimum.\n"
3156 "Be aware that negative eigenvalues may occur\n"
3157 "when the resulting matrix is diagonalized.");
3160 /***********************************************************
3162 * Loop over all pairs in matrix
3164 * do_force called twice. Once with positive and
3165 * once with negative displacement
3167 ************************************************************/
3169 /* Steps are divided one by one over the nodes */
3171 auto state_work_x = makeArrayRef(state_work.s.x);
3172 auto state_work_f = state_work.f.view().force();
3173 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
3175 size_t atom = atom_index[aid];
3176 for (size_t d = 0; d < DIM; d++)
3179 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
3182 x_min = state_work_x[atom][d];
3184 for (unsigned int dx = 0; (dx < 2); dx++)
3188 state_work_x[atom][d] = x_min - der_range;
3192 state_work_x[atom][d] = x_min + der_range;
3195 /* Make evaluate_energy do a single node force calculation */
3199 /* Now is the time to relax the shells */
3200 relax_shell_flexcon(fplog,
3203 mdrunOptions.verbose,
3214 state_work.s.natoms,
3215 state_work.s.x.arrayRefWithPadding(),
3216 state_work.s.v.arrayRefWithPadding(),
3218 state_work.s.lambda,
3220 &state_work.f.view(),
3231 DDBalanceRegionHandler(nullptr));
3237 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid * 2 + dx, FALSE);
3240 cr->nnodes = nnodes;
3244 std::copy(state_work_f.begin(), state_work_f.begin() + atom_index.size(), fneg.begin());
3248 /* x is restored to original */
3249 state_work_x[atom][d] = x_min;
3251 for (size_t j = 0; j < atom_index.size(); j++)
3253 for (size_t k = 0; (k < DIM); k++)
3255 dfdx[j][k] = -(state_work_f[atom_index[j]][k] - fneg[j][k]) / (2 * der_range);
3262 # define mpi_type GMX_MPI_REAL
3263 MPI_Send(dfdx[0], atom_index.size() * DIM, mpi_type, MASTER(cr), cr->nodeid, cr->mpi_comm_mygroup);
3268 for (index node = 0; (node < nnodes && aid + node < ssize(atom_index)); node++)
3274 MPI_Recv(dfdx[0], atom_index.size() * DIM, mpi_type, node, node, cr->mpi_comm_mygroup, &stat);
3279 row = (aid + node) * DIM + d;
3281 for (size_t j = 0; j < atom_index.size(); j++)
3283 for (size_t k = 0; k < DIM; k++)
3289 if (col >= row && dfdx[j][k] != 0.0)
3291 gmx_sparsematrix_increment_value(sparse_matrix, row, col, dfdx[j][k]);
3296 full_matrix[row * sz + col] = dfdx[j][k];
3303 if (mdrunOptions.verbose && fplog)
3308 /* write progress */
3309 if (bIsMaster && mdrunOptions.verbose)
3312 "\rFinished step %d out of %td",
3313 std::min<int>(atom + nnodes, atom_index.size()),
3321 fprintf(stderr, "\n\nWriting Hessian...\n");
3322 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
3325 finish_em(cr, outf, walltime_accounting, wcycle);
3327 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size() * 2);