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>();
397 initialize_lambdas(fplog, *ir, MASTER(cr), fep_state, lambda);
399 if (ir->eI == IntegrationAlgorithm::NM)
401 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
403 *shellfc = init_shell_flexcon(stdout,
405 constr ? constr->numFlexibleConstraints() : 0,
408 thisRankHasDuty(cr, DUTY_PME));
412 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI),
413 "This else currently only handles energy minimizers, consider if your algorithm "
414 "needs shell/flexible-constraint support");
416 /* With energy minimization, shells and flexible constraints are
417 * automatically minimized when treated like normal DOFS.
419 if (shellfc != nullptr)
425 if (DOMAINDECOMP(cr))
427 dd_init_local_state(*cr->dd, state_global, &ems->s);
429 /* Distribute the charge groups over the nodes from the master node */
430 dd_partition_system(fplog,
451 dd_store_state(*cr->dd, &ems->s);
455 state_change_natoms(state_global, state_global->natoms);
456 /* Just copy the state */
457 ems->s = *state_global;
458 state_change_natoms(&ems->s, ems->s.natoms);
460 mdAlgorithmsSetupAtomData(
461 cr, *ir, top_global, top, fr, &ems->f, mdAtoms, constr, vsite, shellfc ? *shellfc : nullptr);
464 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[FreeEnergyPerturbationCouplingType::Mass]);
468 // TODO how should this cross-module support dependency be managed?
469 if (ir->eConstrAlg == ConstraintAlgorithm::Shake && gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
472 "Can not do energy minimization with %s, use %s\n",
473 enumValueToString(ConstraintAlgorithm::Shake),
474 enumValueToString(ConstraintAlgorithm::Lincs));
477 if (!ir->bContinuation)
479 /* Constrain the starting coordinates */
480 bool needsLogging = true;
481 bool computeEnergy = true;
482 bool computeVirial = false;
484 constr->apply(needsLogging,
489 ems->s.x.arrayRefWithPadding(),
490 ems->s.x.arrayRefWithPadding(),
493 ems->s.lambda[FreeEnergyPerturbationCouplingType::Fep],
495 gmx::ArrayRefWithPadding<RVec>(),
498 gmx::ConstraintVariable::Positions);
504 *gstat = global_stat_init(ir);
511 calc_shifts(ems->s.box, fr->shift_vec);
514 //! Finalize the minimization
515 static void finish_em(const t_commrec* cr,
517 gmx_walltime_accounting_t walltime_accounting,
518 gmx_wallcycle_t wcycle)
520 if (!thisRankHasDuty(cr, DUTY_PME))
522 /* Tell the PME only node to finish */
523 gmx_pme_send_finish(cr);
528 em_time_end(walltime_accounting, wcycle);
531 //! Swap two different EM states during minimization
532 static void swap_em_state(em_state_t** ems1, em_state_t** ems2)
541 //! Save the EM trajectory
542 static void write_em_traj(FILE* fplog,
548 const gmx_mtop_t& top_global,
549 const t_inputrec* ir,
552 t_state* state_global,
553 ObservablesHistory* observablesHistory)
559 mdof_flags |= MDOF_X;
563 mdof_flags |= MDOF_F;
566 /* If we want IMD output, set appropriate MDOF flag */
569 mdof_flags |= MDOF_IMD;
572 gmx::WriteCheckpointDataHolder checkpointDataHolder;
573 mdoutf_write_to_trajectory_files(fplog,
579 static_cast<double>(step),
583 state->f.view().force(),
584 &checkpointDataHolder);
586 if (confout != nullptr)
588 if (DOMAINDECOMP(cr))
590 /* If bX=true, x was collected to state_global in the call above */
593 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
595 cr->dd, state->s.ddp_count, state->s.ddp_count_cg_gl, state->s.cg_gl, state->s.x, globalXRef);
600 /* Copy the local state pointer */
601 state_global = &state->s;
606 if (ir->pbcType != PbcType::No && !ir->bPeriodicMols && DOMAINDECOMP(cr))
608 /* Make molecules whole only for confout writing */
609 do_pbc_mtop(ir->pbcType, state->s.box, &top_global, state_global->x.rvec_array());
612 write_sto_conf_mtop(confout,
615 state_global->x.rvec_array(),
623 //! \brief Do one minimization step
625 // \returns true when the step succeeded, false when a constraint error occurred
626 static bool do_em_step(const t_commrec* cr,
627 const t_inputrec* ir,
631 gmx::ArrayRefWithPadding<const gmx::RVec> force,
633 gmx::Constraints* constr,
640 int nthreads gmx_unused;
642 bool validStep = true;
647 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
649 gmx_incons("state mismatch in do_em_step");
652 s2->flags = s1->flags;
654 if (s2->natoms != s1->natoms)
656 state_change_natoms(s2, s1->natoms);
657 ems2->f.resize(s2->natoms);
659 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
661 s2->cg_gl.resize(s1->cg_gl.size());
664 copy_mat(s1->box, s2->box);
665 /* Copy free energy state */
666 s2->lambda = s1->lambda;
667 copy_mat(s1->box, s2->box);
672 nthreads = gmx_omp_nthreads_get(emntUpdate);
673 #pragma omp parallel num_threads(nthreads)
675 const rvec* x1 = s1->x.rvec_array();
676 rvec* x2 = s2->x.rvec_array();
677 const rvec* f = as_rvec_array(force.unpaddedArrayRef().data());
680 #pragma omp for schedule(static) nowait
681 for (int i = start; i < end; i++)
689 for (int m = 0; m < DIM; m++)
691 if (ir->opts.nFreeze[gf][m])
697 x2[i][m] = x1[i][m] + a * f[i][m];
701 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
704 if (s2->flags & enumValueToBitMask(StateEntry::Cgp))
706 /* Copy the CG p vector */
707 const rvec* p1 = s1->cg_p.rvec_array();
708 rvec* p2 = s2->cg_p.rvec_array();
709 #pragma omp for schedule(static) nowait
710 for (int i = start; i < end; i++)
712 // Trivial OpenMP block that does not throw
713 copy_rvec(p1[i], p2[i]);
717 if (DOMAINDECOMP(cr))
719 /* OpenMP does not supported unsigned loop variables */
720 #pragma omp for schedule(static) nowait
721 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
723 s2->cg_gl[i] = s1->cg_gl[i];
728 if (DOMAINDECOMP(cr))
730 s2->ddp_count = s1->ddp_count;
731 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
737 validStep = constr->apply(TRUE,
742 s1->x.arrayRefWithPadding(),
743 s2->x.arrayRefWithPadding(),
746 s2->lambda[FreeEnergyPerturbationCouplingType::Bonded],
748 gmx::ArrayRefWithPadding<RVec>(),
751 gmx::ConstraintVariable::Positions);
755 /* This global reduction will affect performance at high
756 * parallelization, but we can not really avoid it.
757 * But usually EM is not run at high parallelization.
759 int reductionBuffer = static_cast<int>(!validStep);
760 gmx_sumi(1, &reductionBuffer, cr);
761 validStep = (reductionBuffer == 0);
764 // We should move this check to the different minimizers
765 if (!validStep && ir->eI != IntegrationAlgorithm::Steep)
768 "The coordinates could not be constrained. Minimizer '%s' can not handle "
769 "constraint failures, use minimizer '%s' before using '%s'.",
770 enumValueToString(ir->eI),
771 enumValueToString(IntegrationAlgorithm::Steep),
772 enumValueToString(ir->eI));
779 //! Prepare EM for using domain decomposition parallellization
780 static void em_dd_partition_system(FILE* fplog,
781 const gmx::MDLogger& mdlog,
784 const gmx_mtop_t& top_global,
785 const t_inputrec* ir,
786 gmx::ImdSession* imdSession,
790 gmx::MDAtoms* mdAtoms,
792 VirtualSitesHandler* vsite,
793 gmx::Constraints* constr,
795 gmx_wallcycle_t wcycle)
797 /* Repartition the domain decomposition */
798 dd_partition_system(fplog,
819 dd_store_state(*cr->dd, &ems->s);
825 /*! \brief Class to handle the work of setting and doing an energy evaluation.
827 * This class is a mere aggregate of parameters to pass to evaluate an
828 * energy, so that future changes to names and types of them consume
829 * less time when refactoring other code.
831 * Aggregate initialization is used, for which the chief risk is that
832 * if a member is added at the end and not all initializer lists are
833 * updated, then the member will be value initialized, which will
834 * typically mean initialization to zero.
836 * Use a braced initializer list to construct one of these. */
837 class EnergyEvaluator
840 /*! \brief Evaluates an energy on the state in \c ems.
842 * \todo In practice, the same objects mu_tot, vir, and pres
843 * are always passed to this function, so we would rather have
844 * them as data members. However, their C-array types are
845 * unsuited for aggregate initialization. When the types
846 * improve, the call signature of this method can be reduced.
848 void run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst);
849 //! Handles logging (deprecated).
852 const gmx::MDLogger& mdlog;
853 //! Handles communication.
855 //! Coordinates multi-simulations.
856 const gmx_multisim_t* ms;
857 //! Holds the simulation topology.
858 const gmx_mtop_t& top_global;
859 //! Holds the domain topology.
861 //! User input options.
862 const t_inputrec* inputrec;
863 //! The Interactive Molecular Dynamics session.
864 gmx::ImdSession* imdSession;
865 //! The pull work object.
867 //! Manages flop accounting.
869 //! Manages wall cycle accounting.
870 gmx_wallcycle_t wcycle;
871 //! Coordinates global reduction.
872 gmx_global_stat_t gstat;
873 //! Handles virtual sites.
874 VirtualSitesHandler* vsite;
875 //! Handles constraints.
876 gmx::Constraints* constr;
877 //! Per-atom data for this domain.
878 gmx::MDAtoms* mdAtoms;
879 //! Handles how to calculate the forces.
881 //! Schedule of force-calculation work each step for this task.
882 MdrunScheduleWorkload* runScheduleWork;
883 //! Stores the computed energies.
884 gmx_enerdata_t* enerd;
887 void EnergyEvaluator::run(em_state_t* ems, rvec mu_tot, tensor vir, tensor pres, int64_t count, gmx_bool bFirst)
891 tensor force_vir, shake_vir, ekin;
895 /* Set the time to the initial time, the time does not change during EM */
896 t = inputrec->init_t;
898 if (bFirst || (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
900 /* This is the first state or an old state used before the last ns */
906 if (inputrec->nstlist > 0)
914 vsite->construct(ems->s.x, {}, ems->s.box, gmx::VSiteOperation::Positions);
917 if (DOMAINDECOMP(cr) && bNS)
919 /* Repartition the domain decomposition */
920 em_dd_partition_system(
921 fplog, mdlog, count, cr, top_global, inputrec, imdSession, pull_work, ems, top, mdAtoms, fr, vsite, constr, nrnb, wcycle);
924 /* Calc force & energy on new trial position */
925 /* do_force always puts the charge groups in the box and shifts again
926 * We do not unshift, so molecules are always whole in congrad.c
941 ems->s.x.arrayRefWithPadding(),
954 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY
955 | (bNS ? GMX_FORCE_NS : 0),
956 DDBalanceRegionHandler(cr));
958 /* Clear the unused shake virial and pressure */
959 clear_mat(shake_vir);
962 /* Communicate stuff when parallel */
963 if (PAR(cr) && inputrec->eI != IntegrationAlgorithm::NM)
965 wallcycle_start(wcycle, ewcMoveE);
974 gmx::ArrayRef<real>{},
976 std::vector<real>(1, terminate),
978 CGLO_ENERGY | CGLO_PRESSURE | CGLO_CONSTRAINT);
980 wallcycle_stop(wcycle, ewcMoveE);
983 if (fr->dispersionCorrection)
985 /* Calculate long range corrections to pressure and energy */
986 const DispersionCorrection::Correction correction = fr->dispersionCorrection->calculate(
987 ems->s.box, ems->s.lambda[FreeEnergyPerturbationCouplingType::Vdw]);
989 enerd->term[F_DISPCORR] = correction.energy;
990 enerd->term[F_EPOT] += correction.energy;
991 enerd->term[F_PRES] += correction.pressure;
992 enerd->term[F_DVDL] += correction.dvdl;
996 enerd->term[F_DISPCORR] = 0;
999 ems->epot = enerd->term[F_EPOT];
1003 /* Project out the constraint components of the force */
1004 bool needsLogging = false;
1005 bool computeEnergy = false;
1006 bool computeVirial = true;
1008 auto f = ems->f.view().forceWithPadding();
1009 constr->apply(needsLogging,
1014 ems->s.x.arrayRefWithPadding(),
1016 f.unpaddedArrayRef(),
1018 ems->s.lambda[FreeEnergyPerturbationCouplingType::Bonded],
1020 gmx::ArrayRefWithPadding<RVec>(),
1023 gmx::ConstraintVariable::ForceDispl);
1024 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
1025 m_add(force_vir, shake_vir, vir);
1029 copy_mat(force_vir, vir);
1033 enerd->term[F_PRES] = calc_pres(fr->pbcType, inputrec->nwall, ems->s.box, ekin, vir, pres);
1035 if (inputrec->efep != FreeEnergyPerturbationType::No)
1037 accumulateKineticLambdaComponents(enerd, ems->s.lambda, *inputrec->fepvals);
1040 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
1042 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
1048 //! Parallel utility summing energies and forces
1049 static double reorder_partsum(const t_commrec* cr,
1050 const t_grpopts* opts,
1051 const gmx_mtop_t& top_global,
1052 const em_state_t* s_min,
1053 const em_state_t* s_b)
1057 fprintf(debug, "Doing reorder_partsum\n");
1060 auto fm = s_min->f.view().force();
1061 auto fb = s_b->f.view().force();
1063 /* Collect fm in a global vector fmg.
1064 * This conflicts with the spirit of domain decomposition,
1065 * but to fully optimize this a much more complicated algorithm is required.
1067 const int natoms = top_global.natoms;
1071 gmx::ArrayRef<const int> indicesMin = s_min->s.cg_gl;
1073 for (int a : indicesMin)
1075 copy_rvec(fm[i], fmg[a]);
1078 gmx_sum(top_global.natoms * 3, fmg[0], cr);
1080 /* Now we will determine the part of the sum for the cgs in state s_b */
1081 gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
1086 gmx::ArrayRef<const unsigned char> grpnrFREEZE =
1087 top_global.groups.groupNumbers[SimulationAtomGroupType::Freeze];
1088 for (int a : indicesB)
1090 if (!grpnrFREEZE.empty())
1092 gf = grpnrFREEZE[i];
1094 for (int m = 0; m < DIM; m++)
1096 if (!opts->nFreeze[gf][m])
1098 partsum += (fb[i][m] - fmg[a][m]) * fb[i][m];
1109 //! Print some stuff, like beta, whatever that means.
1110 static real pr_beta(const t_commrec* cr,
1111 const t_grpopts* opts,
1113 const gmx_mtop_t& top_global,
1114 const em_state_t* s_min,
1115 const em_state_t* s_b)
1119 /* This is just the classical Polak-Ribiere calculation of beta;
1120 * it looks a bit complicated since we take freeze groups into account,
1121 * and might have to sum it in parallel runs.
1124 if (!DOMAINDECOMP(cr)
1125 || (s_min->s.ddp_count == cr->dd->ddp_count && s_b->s.ddp_count == cr->dd->ddp_count))
1127 auto fm = s_min->f.view().force();
1128 auto fb = s_b->f.view().force();
1131 /* This part of code can be incorrect with DD,
1132 * since the atom ordering in s_b and s_min might differ.
1134 for (int i = 0; i < mdatoms->homenr; i++)
1136 if (mdatoms->cFREEZE)
1138 gf = mdatoms->cFREEZE[i];
1140 for (int m = 0; m < DIM; m++)
1142 if (!opts->nFreeze[gf][m])
1144 sum += (fb[i][m] - fm[i][m]) * fb[i][m];
1151 /* We need to reorder cgs while summing */
1152 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1156 gmx_sumd(1, &sum, cr);
1159 return sum / gmx::square(s_min->fnorm);
1165 void LegacySimulator::do_cg()
1167 const char* CG = "Polak-Ribiere Conjugate Gradients";
1169 gmx_localtop_t top(top_global.ffparams);
1170 gmx_global_stat_t gstat;
1171 double tmp, minstep;
1173 real a, b, c, beta = 0.0;
1176 gmx_bool converged, foundlower;
1177 rvec mu_tot = { 0 };
1178 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1180 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1181 int m, step, nminstep;
1182 auto mdatoms = mdAtoms->mdatoms();
1187 "Note that activating conjugate gradient energy minimization via the "
1188 "integrator .mdp option and the command gmx mdrun may "
1189 "be available in a different form in a future version of GROMACS, "
1190 "e.g. gmx minimize and an .mdp option.");
1196 // In CG, the state is extended with a search direction
1197 state_global->flags |= enumValueToBitMask(StateEntry::Cgp);
1199 // Ensure the extra per-atom state array gets allocated
1200 state_change_natoms(state_global, state_global->natoms);
1202 // Initialize the search direction to zero
1203 for (RVec& cg_p : state_global->cg_p)
1209 /* Create 4 states on the stack and extract pointers that we will swap */
1210 em_state_t s0{}, s1{}, s2{}, s3{};
1211 em_state_t* s_min = &s0;
1212 em_state_t* s_a = &s1;
1213 em_state_t* s_b = &s2;
1214 em_state_t* s_c = &s3;
1216 /* Init em and store the local state in s_min */
1235 const bool simulationsShareState = false;
1236 gmx_mdoutf* outf = init_mdoutf(fplog,
1247 StartingBehavior::NewSimulation,
1248 simulationsShareState,
1250 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
1256 StartingBehavior::NewSimulation,
1257 simulationsShareState,
1260 /* Print to log file */
1261 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1263 /* Max number of steps */
1264 number_steps = inputrec->nsteps;
1268 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1272 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1275 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
1276 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1277 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
1278 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1279 /* do_force always puts the charge groups in the box and shifts again
1280 * We do not unshift, so molecules are always whole in congrad.c
1282 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1286 /* Copy stuff to the energy bin for easy printing etc. */
1287 matrix nullBox = {};
1288 energyOutput.addDataAtEnergyStep(false,
1290 static_cast<double>(step),
1306 EnergyOutput::printHeader(fplog, step, step);
1307 energyOutput.printStepToEnergyFile(
1308 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
1311 /* Estimate/guess the initial stepsize */
1312 stepsize = inputrec->em_stepsize / s_min->fnorm;
1316 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1317 fprintf(stderr, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1318 fprintf(stderr, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1319 fprintf(stderr, "\n");
1320 /* and copy to the log file too... */
1321 fprintf(fplog, " F-max = %12.5e on atom %d\n", s_min->fmax, s_min->a_fmax + 1);
1322 fprintf(fplog, " F-Norm = %12.5e\n", s_min->fnorm / sqrtNumAtoms);
1323 fprintf(fplog, "\n");
1325 /* Start the loop over CG steps.
1326 * Each successful step is counted, and we continue until
1327 * we either converge or reach the max number of steps.
1330 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1333 /* start taking steps in a new direction
1334 * First time we enter the routine, beta=0, and the direction is
1335 * simply the negative gradient.
1338 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1339 gmx::ArrayRef<gmx::RVec> pm = s_min->s.cg_p;
1340 gmx::ArrayRef<const gmx::RVec> sfm = s_min->f.view().force();
1343 for (int i = 0; i < mdatoms->homenr; i++)
1345 if (mdatoms->cFREEZE)
1347 gf = mdatoms->cFREEZE[i];
1349 for (m = 0; m < DIM; m++)
1351 if (!inputrec->opts.nFreeze[gf][m])
1353 pm[i][m] = sfm[i][m] + beta * pm[i][m];
1354 gpa -= pm[i][m] * sfm[i][m];
1355 /* f is negative gradient, thus the sign */
1364 /* Sum the gradient along the line across CPUs */
1367 gmx_sumd(1, &gpa, cr);
1370 /* Calculate the norm of the search vector */
1371 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1373 /* Just in case stepsize reaches zero due to numerical precision... */
1376 stepsize = inputrec->em_stepsize / pnorm;
1380 * Double check the value of the derivative in the search direction.
1381 * If it is positive it must be due to the old information in the
1382 * CG formula, so just remove that and start over with beta=0.
1383 * This corresponds to a steepest descent step.
1388 step--; /* Don't count this step since we are restarting */
1389 continue; /* Go back to the beginning of the big for-loop */
1392 /* Calculate minimum allowed stepsize, before the average (norm)
1393 * relative change in coordinate is smaller than precision
1396 auto s_min_x = makeArrayRef(s_min->s.x);
1397 for (int i = 0; i < mdatoms->homenr; i++)
1399 for (m = 0; m < DIM; m++)
1401 tmp = fabs(s_min_x[i][m]);
1406 tmp = pm[i][m] / tmp;
1407 minstep += tmp * tmp;
1410 /* Add up from all CPUs */
1413 gmx_sumd(1, &minstep, cr);
1416 minstep = GMX_REAL_EPS / sqrt(minstep / (3 * top_global.natoms));
1418 if (stepsize < minstep)
1424 /* Write coordinates if necessary */
1425 do_x = do_per_step(step, inputrec->nstxout);
1426 do_f = do_per_step(step, inputrec->nstfout);
1429 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, step, s_min, state_global, observablesHistory);
1431 /* Take a step downhill.
1432 * In theory, we should minimize the function along this direction.
1433 * That is quite possible, but it turns out to take 5-10 function evaluations
1434 * for each line. However, we dont really need to find the exact minimum -
1435 * it is much better to start a new CG step in a modified direction as soon
1436 * as we are close to it. This will save a lot of energy evaluations.
1438 * In practice, we just try to take a single step.
1439 * If it worked (i.e. lowered the energy), we increase the stepsize but
1440 * the continue straight to the next CG step without trying to find any minimum.
1441 * If it didn't work (higher energy), there must be a minimum somewhere between
1442 * the old position and the new one.
1444 * Due to the finite numerical accuracy, it turns out that it is a good idea
1445 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1446 * This leads to lower final energies in the tests I've done. / Erik
1448 s_a->epot = s_min->epot;
1450 c = a + stepsize; /* reference position along line is zero */
1452 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1454 em_dd_partition_system(fplog,
1472 /* Take a trial step (new coords in s_c) */
1473 do_em_step(cr, inputrec, mdatoms, s_min, c, s_min->s.cg_p.constArrayRefWithPadding(), s_c, constr, -1);
1476 /* Calculate energy for the trial step */
1477 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1479 /* Calc derivative along line */
1480 const rvec* pc = s_c->s.cg_p.rvec_array();
1481 gmx::ArrayRef<const gmx::RVec> sfc = s_c->f.view().force();
1483 for (int i = 0; i < mdatoms->homenr; i++)
1485 for (m = 0; m < DIM; m++)
1487 gpc -= pc[i][m] * sfc[i][m]; /* f is negative gradient, thus the sign */
1490 /* Sum the gradient along the line across CPUs */
1493 gmx_sumd(1, &gpc, cr);
1496 /* This is the max amount of increase in energy we tolerate */
1497 tmp = std::sqrt(GMX_REAL_EPS) * fabs(s_a->epot);
1499 /* Accept the step if the energy is lower, or if it is not significantly higher
1500 * and the line derivative is still negative.
1502 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1505 /* Great, we found a better energy. Increase step for next iteration
1506 * if we are still going down, decrease it otherwise
1510 stepsize *= 1.618034; /* The golden section */
1514 stepsize *= 0.618034; /* 1/golden section */
1519 /* New energy is the same or higher. We will have to do some work
1520 * to find a smaller value in the interval. Take smaller step next time!
1523 stepsize *= 0.618034;
1527 /* OK, if we didn't find a lower value we will have to locate one now - there must
1528 * be one in the interval [a=0,c].
1529 * The same thing is valid here, though: Don't spend dozens of iterations to find
1530 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1531 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1533 * I also have a safeguard for potentially really pathological functions so we never
1534 * take more than 20 steps before we give up ...
1536 * If we already found a lower value we just skip this step and continue to the update.
1545 /* Select a new trial point.
1546 * If the derivatives at points a & c have different sign we interpolate to zero,
1547 * otherwise just do a bisection.
1549 if (gpa < 0 && gpc > 0)
1551 b = a + gpa * (a - c) / (gpc - gpa);
1558 /* safeguard if interpolation close to machine accuracy causes errors:
1559 * never go outside the interval
1561 if (b <= a || b >= c)
1566 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1568 /* Reload the old state */
1569 em_dd_partition_system(fplog,
1587 /* Take a trial step to this new point - new coords in s_b */
1588 do_em_step(cr, inputrec, mdatoms, s_min, b, s_min->s.cg_p.constArrayRefWithPadding(), s_b, constr, -1);
1591 /* Calculate energy for the trial step */
1592 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1594 /* p does not change within a step, but since the domain decomposition
1595 * might change, we have to use cg_p of s_b here.
1597 const rvec* pb = s_b->s.cg_p.rvec_array();
1598 gmx::ArrayRef<const gmx::RVec> sfb = s_b->f.view().force();
1600 for (int i = 0; i < mdatoms->homenr; i++)
1602 for (m = 0; m < DIM; m++)
1604 gpb -= pb[i][m] * sfb[i][m]; /* f is negative gradient, thus the sign */
1607 /* Sum the gradient along the line across CPUs */
1610 gmx_sumd(1, &gpb, cr);
1615 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a->epot, s_b->epot, s_c->epot, gpb);
1618 epot_repl = s_b->epot;
1620 /* Keep one of the intervals based on the value of the derivative at the new point */
1623 /* Replace c endpoint with b */
1624 swap_em_state(&s_b, &s_c);
1630 /* Replace a endpoint with b */
1631 swap_em_state(&s_b, &s_a);
1637 * Stop search as soon as we find a value smaller than the endpoints.
1638 * Never run more than 20 steps, no matter what.
1641 } while ((epot_repl > s_a->epot || epot_repl > s_c->epot) && (nminstep < 20));
1643 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot) * GMX_REAL_EPS || nminstep >= 20)
1645 /* OK. We couldn't find a significantly lower energy.
1646 * If beta==0 this was steepest descent, and then we give up.
1647 * If not, set beta=0 and restart with steepest descent before quitting.
1657 /* Reset memory before giving up */
1663 /* Select min energy state of A & C, put the best in B.
1665 if (s_c->epot < s_a->epot)
1669 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c->epot, s_a->epot);
1671 swap_em_state(&s_b, &s_c);
1678 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a->epot, s_c->epot);
1680 swap_em_state(&s_b, &s_a);
1688 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n", s_c->epot);
1690 swap_em_state(&s_b, &s_c);
1694 /* new search direction */
1695 /* beta = 0 means forget all memory and restart with steepest descents. */
1696 if (nstcg && ((step % nstcg) == 0))
1702 /* s_min->fnorm cannot be zero, because then we would have converged
1706 /* Polak-Ribiere update.
1707 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1709 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1711 /* Limit beta to prevent oscillations */
1712 if (fabs(beta) > 5.0)
1718 /* update positions */
1719 swap_em_state(&s_min, &s_b);
1722 /* Print it if necessary */
1725 if (mdrunOptions.verbose)
1727 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1729 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1732 s_min->fnorm / sqrtNumAtoms,
1737 /* Store the new (lower) energies */
1738 matrix nullBox = {};
1739 energyOutput.addDataAtEnergyStep(false,
1741 static_cast<double>(step),
1757 do_log = do_per_step(step, inputrec->nstlog);
1758 do_ene = do_per_step(step, inputrec->nstenergy);
1760 imdSession->fillEnergyRecord(step, TRUE);
1764 EnergyOutput::printHeader(fplog, step, step);
1766 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1770 do_log ? fplog : nullptr,
1777 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1778 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x, 0))
1780 imdSession->sendPositionsAndEnergies();
1783 /* Stop when the maximum force lies below tolerance.
1784 * If we have reached machine precision, converged is already set to true.
1786 converged = converged || (s_min->fmax < inputrec->em_tol);
1788 } /* End of the loop */
1792 step--; /* we never took that last step in this case */
1794 if (s_min->fmax > inputrec->em_tol)
1798 warn_step(fplog, inputrec->em_tol, s_min->fmax, step - 1 == number_steps, FALSE);
1805 /* If we printed energy and/or logfile last step (which was the last step)
1806 * we don't have to do it again, but otherwise print the final values.
1810 /* Write final value to log since we didn't do anything the last step */
1811 EnergyOutput::printHeader(fplog, step, step);
1813 if (!do_ene || !do_log)
1815 /* Write final energy file entries */
1816 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
1820 !do_log ? fplog : nullptr,
1828 /* Print some stuff... */
1831 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1835 * For accurate normal mode calculation it is imperative that we
1836 * store the last conformation into the full precision binary trajectory.
1838 * However, we should only do it if we did NOT already write this step
1839 * above (which we did if do_x or do_f was true).
1841 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1842 * in the trajectory with the same step number.
1844 do_x = !do_per_step(step, inputrec->nstxout);
1845 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1848 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, s_min, state_global, observablesHistory);
1853 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1854 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1855 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps, s_min, sqrtNumAtoms);
1857 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1860 finish_em(cr, outf, walltime_accounting, wcycle);
1862 /* To print the actual number of steps we needed somewhere */
1863 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1867 void LegacySimulator::do_lbfgs()
1869 static const char* LBFGS = "Low-Memory BFGS Minimizer";
1871 gmx_localtop_t top(top_global.ffparams);
1872 gmx_global_stat_t gstat;
1873 auto mdatoms = mdAtoms->mdatoms();
1878 "Note that activating L-BFGS energy minimization via the "
1879 "integrator .mdp option and the command gmx mdrun may "
1880 "be available in a different form in a future version of GROMACS, "
1881 "e.g. gmx minimize and an .mdp option.");
1885 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1888 if (nullptr != constr)
1892 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1893 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1896 const int n = 3 * state_global->natoms;
1897 const int nmaxcorr = inputrec->nbfgscorr;
1899 std::vector<real> p(n);
1900 std::vector<real> rho(nmaxcorr);
1901 std::vector<real> alpha(nmaxcorr);
1903 std::vector<std::vector<real>> dx(nmaxcorr);
1904 for (auto& dxCorr : dx)
1909 std::vector<std::vector<real>> dg(nmaxcorr);
1910 for (auto& dgCorr : dg)
1937 const bool simulationsShareState = false;
1938 gmx_mdoutf* outf = init_mdoutf(fplog,
1949 StartingBehavior::NewSimulation,
1950 simulationsShareState,
1952 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
1958 StartingBehavior::NewSimulation,
1959 simulationsShareState,
1962 const int start = 0;
1963 const int end = mdatoms->homenr;
1965 /* We need 4 working states */
1966 em_state_t s0{}, s1{}, s2{}, s3{};
1967 em_state_t* sa = &s0;
1968 em_state_t* sb = &s1;
1969 em_state_t* sc = &s2;
1970 em_state_t* last = &s3;
1971 /* Initialize by copying the state from ems (we could skip x and f here) */
1976 /* Print to log file */
1977 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1979 /* Max number of steps */
1980 const int number_steps = inputrec->nsteps;
1982 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1983 std::vector<bool> frozen(n);
1985 for (int i = start; i < end; i++)
1987 if (mdatoms->cFREEZE)
1989 gf = mdatoms->cFREEZE[i];
1991 for (int m = 0; m < DIM; m++)
1993 frozen[3 * i + m] = (inputrec->opts.nFreeze[gf][m] != 0);
1998 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
2002 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
2007 vsite->construct(state_global->x, {}, state_global->box, VSiteOperation::Positions);
2010 /* Call the force routine and some auxiliary (neighboursearching etc.) */
2011 /* do_force always puts the charge groups in the box and shifts again
2012 * We do not unshift, so molecules are always whole
2015 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2016 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2017 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2021 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
2025 /* Copy stuff to the energy bin for easy printing etc. */
2026 matrix nullBox = {};
2027 energyOutput.addDataAtEnergyStep(false,
2029 static_cast<double>(step),
2045 EnergyOutput::printHeader(fplog, step, step);
2046 energyOutput.printStepToEnergyFile(
2047 mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, fr->fcdata.get(), nullptr);
2050 /* Set the initial step.
2051 * since it will be multiplied by the non-normalized search direction
2052 * vector (force vector the first time), we scale it by the
2053 * norm of the force.
2058 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2059 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2060 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2061 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2062 fprintf(stderr, "\n");
2063 /* and copy to the log file too... */
2064 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
2065 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
2066 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm / sqrtNumAtoms);
2067 fprintf(fplog, "\n");
2070 // Point is an index to the memory of search directions, where 0 is the first one.
2073 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
2074 real* fInit = static_cast<real*>(ems.f.view().force().data()[0]);
2075 for (int i = 0; i < n; i++)
2079 dx[point][i] = fInit[i]; /* Initial search direction */
2087 // Stepsize will be modified during the search, and actually it is not critical
2088 // (the main efficiency in the algorithm comes from changing directions), but
2089 // we still need an initial value, so estimate it as the inverse of the norm
2090 // so we take small steps where the potential fluctuates a lot.
2091 double stepsize = 1.0 / ems.fnorm;
2093 /* Start the loop over BFGS steps.
2094 * Each successful step is counted, and we continue until
2095 * we either converge or reach the max number of steps.
2103 /* Set the gradient from the force */
2104 bool converged = false;
2105 for (int step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
2108 /* Write coordinates if necessary */
2109 const bool do_x = do_per_step(step, inputrec->nstxout);
2110 const bool do_f = do_per_step(step, inputrec->nstfout);
2115 mdof_flags |= MDOF_X;
2120 mdof_flags |= MDOF_F;
2125 mdof_flags |= MDOF_IMD;
2128 gmx::WriteCheckpointDataHolder checkpointDataHolder;
2129 mdoutf_write_to_trajectory_files(fplog,
2135 static_cast<real>(step),
2139 ems.f.view().force(),
2140 &checkpointDataHolder);
2142 /* Do the linesearching in the direction dx[point][0..(n-1)] */
2144 /* make s a pointer to current search direction - point=0 first time we get here */
2145 gmx::ArrayRef<const real> s = dx[point];
2147 const real* xx = static_cast<real*>(ems.s.x.rvec_array()[0]);
2148 const real* ff = static_cast<real*>(ems.f.view().force().data()[0]);
2150 // calculate line gradient in position A
2152 for (int i = 0; i < n; i++)
2154 gpa -= s[i] * ff[i];
2157 /* Calculate minimum allowed stepsize along the line, before the average (norm)
2158 * relative change in coordinate is smaller than precision
2161 for (int i = 0; i < n; i++)
2163 double tmp = fabs(xx[i]);
2169 minstep += tmp * tmp;
2171 minstep = GMX_REAL_EPS / sqrt(minstep / n);
2173 if (stepsize < minstep)
2179 // Before taking any steps along the line, store the old position
2181 real* lastx = static_cast<real*>(last->s.x.data()[0]);
2182 real* lastf = static_cast<real*>(last->f.view().force().data()[0]);
2183 const real Epot0 = ems.epot;
2187 /* Take a step downhill.
2188 * In theory, we should find the actual minimum of the function in this
2189 * direction, somewhere along the line.
2190 * That is quite possible, but it turns out to take 5-10 function evaluations
2191 * for each line. However, we dont really need to find the exact minimum -
2192 * it is much better to start a new BFGS step in a modified direction as soon
2193 * as we are close to it. This will save a lot of energy evaluations.
2195 * In practice, we just try to take a single step.
2196 * If it worked (i.e. lowered the energy), we increase the stepsize but
2197 * continue straight to the next BFGS step without trying to find any minimum,
2198 * i.e. we change the search direction too. If the line was smooth, it is
2199 * likely we are in a smooth region, and then it makes sense to take longer
2200 * steps in the modified search direction too.
2202 * If it didn't work (higher energy), there must be a minimum somewhere between
2203 * the old position and the new one. Then we need to start by finding a lower
2204 * value before we change search direction. Since the energy was apparently
2205 * quite rough, we need to decrease the step size.
2207 * Due to the finite numerical accuracy, it turns out that it is a good idea
2208 * to accept a SMALL increase in energy, if the derivative is still downhill.
2209 * This leads to lower final energies in the tests I've done. / Erik
2212 // State "A" is the first position along the line.
2213 // reference position along line is initially zero
2216 // Check stepsize first. We do not allow displacements
2217 // larger than emstep.
2223 // Pick a new position C by adding stepsize to A.
2226 // Calculate what the largest change in any individual coordinate
2227 // would be (translation along line * gradient along line)
2229 for (int i = 0; i < n; i++)
2231 real delta = c * s[i];
2232 if (delta > maxdelta)
2237 // If any displacement is larger than the stepsize limit, reduce the step
2238 if (maxdelta > inputrec->em_stepsize)
2242 } while (maxdelta > inputrec->em_stepsize);
2244 // Take a trial step and move the coordinate array xc[] to position C
2245 real* xc = static_cast<real*>(sc->s.x.rvec_array()[0]);
2246 for (int i = 0; i < n; i++)
2248 xc[i] = lastx[i] + c * s[i];
2252 // Calculate energy for the trial step in position C
2253 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2255 // Calc line gradient in position C
2256 real* fc = static_cast<real*>(sc->f.view().force()[0]);
2258 for (int i = 0; i < n; i++)
2260 gpc -= s[i] * fc[i]; /* f is negative gradient, thus the sign */
2262 /* Sum the gradient along the line across CPUs */
2265 gmx_sumd(1, &gpc, cr);
2268 // This is the max amount of increase in energy we tolerate.
2269 // By allowing VERY small changes (close to numerical precision) we
2270 // frequently find even better (lower) final energies.
2271 double tmp = std::sqrt(GMX_REAL_EPS) * fabs(sa->epot);
2273 // Accept the step if the energy is lower in the new position C (compared to A),
2274 // or if it is not significantly higher and the line derivative is still negative.
2275 bool foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2276 // If true, great, we found a better energy. We no longer try to alter the
2277 // stepsize, but simply accept this new better position. The we select a new
2278 // search direction instead, which will be much more efficient than continuing
2279 // to take smaller steps along a line. Set fnorm based on the new C position,
2280 // which will be used to update the stepsize to 1/fnorm further down.
2282 // If false, the energy is NOT lower in point C, i.e. it will be the same
2283 // or higher than in point A. In this case it is pointless to move to point C,
2284 // so we will have to do more iterations along the same line to find a smaller
2285 // value in the interval [A=0.0,C].
2286 // Here, A is still 0.0, but that will change when we do a search in the interval
2287 // [0.0,C] below. That search we will do by interpolation or bisection rather
2288 // than with the stepsize, so no need to modify it. For the next search direction
2289 // it will be reset to 1/fnorm anyway.
2294 // OK, if we didn't find a lower value we will have to locate one now - there must
2295 // be one in the interval [a,c].
2296 // The same thing is valid here, though: Don't spend dozens of iterations to find
2297 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2298 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2299 // I also have a safeguard for potentially really pathological functions so we never
2300 // take more than 20 steps before we give up.
2301 // If we already found a lower value we just skip this step and continue to the update.
2306 // Select a new trial point B in the interval [A,C].
2307 // If the derivatives at points a & c have different sign we interpolate to zero,
2308 // otherwise just do a bisection since there might be multiple minima/maxima
2309 // inside the interval.
2311 if (gpa < 0 && gpc > 0)
2313 b = a + gpa * (a - c) / (gpc - gpa);
2320 /* safeguard if interpolation close to machine accuracy causes errors:
2321 * never go outside the interval
2323 if (b <= a || b >= c)
2328 // Take a trial step to point B
2329 real* xb = static_cast<real*>(sb->s.x.rvec_array()[0]);
2330 for (int i = 0; i < n; i++)
2332 xb[i] = lastx[i] + b * s[i];
2336 // Calculate energy for the trial step in point B
2337 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2340 // Calculate gradient in point B
2341 real* fb = static_cast<real*>(sb->f.view().force()[0]);
2343 for (int i = 0; i < n; i++)
2345 gpb -= s[i] * fb[i]; /* f is negative gradient, thus the sign */
2347 /* Sum the gradient along the line across CPUs */
2350 gmx_sumd(1, &gpb, cr);
2353 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2354 // at the new point B, and rename the endpoints of this new interval A and C.
2357 /* Replace c endpoint with b */
2359 /* copy state b to c */
2364 /* Replace a endpoint with b */
2366 /* copy state b to a */
2371 * Stop search as soon as we find a value smaller than the endpoints,
2372 * or if the tolerance is below machine precision.
2373 * Never run more than 20 steps, no matter what.
2376 } while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2378 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2380 /* OK. We couldn't find a significantly lower energy.
2381 * If ncorr==0 this was steepest descent, and then we give up.
2382 * If not, reset memory to restart as steepest descent before quitting.
2394 /* Search in gradient direction */
2395 for (int i = 0; i < n; i++)
2397 dx[point][i] = ff[i];
2399 /* Reset stepsize */
2400 stepsize = 1.0 / fnorm;
2405 /* Select min energy state of A & C, put the best in xx/ff/Epot
2407 if (sc->epot < sa->epot)
2428 /* Update the memory information, and calculate a new
2429 * approximation of the inverse hessian
2432 /* Have new data in Epot, xx, ff */
2433 if (ncorr < nmaxcorr)
2438 for (int i = 0; i < n; i++)
2440 dg[point][i] = lastf[i] - ff[i];
2441 dx[point][i] *= step_taken;
2446 for (int i = 0; i < n; i++)
2448 dgdg += dg[point][i] * dg[point][i];
2449 dgdx += dg[point][i] * dx[point][i];
2452 const real diag = dgdx / dgdg;
2454 rho[point] = 1.0 / dgdx;
2457 if (point >= nmaxcorr)
2463 for (int i = 0; i < n; i++)
2470 /* Recursive update. First go back over the memory points */
2471 for (int k = 0; k < ncorr; k++)
2480 for (int i = 0; i < n; i++)
2482 sq += dx[cp][i] * p[i];
2485 alpha[cp] = rho[cp] * sq;
2487 for (int i = 0; i < n; i++)
2489 p[i] -= alpha[cp] * dg[cp][i];
2493 for (int i = 0; i < n; i++)
2498 /* And then go forward again */
2499 for (int k = 0; k < ncorr; k++)
2502 for (int i = 0; i < n; i++)
2504 yr += p[i] * dg[cp][i];
2507 real beta = rho[cp] * yr;
2508 beta = alpha[cp] - beta;
2510 for (int i = 0; i < n; i++)
2512 p[i] += beta * dx[cp][i];
2522 for (int i = 0; i < n; i++)
2526 dx[point][i] = p[i];
2534 /* Print it if necessary */
2537 if (mdrunOptions.verbose)
2539 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2541 "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2544 ems.fnorm / sqrtNumAtoms,
2549 /* Store the new (lower) energies */
2550 matrix nullBox = {};
2551 energyOutput.addDataAtEnergyStep(false,
2553 static_cast<double>(step),
2569 do_log = do_per_step(step, inputrec->nstlog);
2570 do_ene = do_per_step(step, inputrec->nstenergy);
2572 imdSession->fillEnergyRecord(step, TRUE);
2576 EnergyOutput::printHeader(fplog, step, step);
2578 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2582 do_log ? fplog : nullptr,
2589 /* Send x and E to IMD client, if bIMD is TRUE. */
2590 if (imdSession->run(step, TRUE, state_global->box, state_global->x, 0) && MASTER(cr))
2592 imdSession->sendPositionsAndEnergies();
2595 // Reset stepsize in we are doing more iterations
2598 /* Stop when the maximum force lies below tolerance.
2599 * If we have reached machine precision, converged is already set to true.
2601 converged = converged || (ems.fmax < inputrec->em_tol);
2603 } /* End of the loop */
2607 step--; /* we never took that last step in this case */
2609 if (ems.fmax > inputrec->em_tol)
2613 warn_step(fplog, inputrec->em_tol, ems.fmax, step - 1 == number_steps, FALSE);
2618 /* If we printed energy and/or logfile last step (which was the last step)
2619 * we don't have to do it again, but otherwise print the final values.
2621 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2623 EnergyOutput::printHeader(fplog, step, step);
2625 if (!do_ene || !do_log) /* Write final energy file entries */
2627 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf),
2631 !do_log ? fplog : nullptr,
2638 /* Print some stuff... */
2641 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2645 * For accurate normal mode calculation it is imperative that we
2646 * store the last conformation into the full precision binary trajectory.
2648 * However, we should only do it if we did NOT already write this step
2649 * above (which we did if do_x or do_f was true).
2651 const bool do_x = !do_per_step(step, inputrec->nstxout);
2652 const bool do_f = !do_per_step(step, inputrec->nstfout);
2654 fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm), top_global, inputrec, step, &ems, state_global, observablesHistory);
2658 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2659 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2660 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged, number_steps, &ems, sqrtNumAtoms);
2662 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2665 finish_em(cr, outf, walltime_accounting, wcycle);
2667 /* To print the actual number of steps we needed somewhere */
2668 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2671 void LegacySimulator::do_steep()
2673 const char* SD = "Steepest Descents";
2674 gmx_localtop_t top(top_global.ffparams);
2675 gmx_global_stat_t gstat;
2678 gmx_bool bDone, bAbort, do_x, do_f;
2680 rvec mu_tot = { 0 };
2683 int steps_accepted = 0;
2684 auto mdatoms = mdAtoms->mdatoms();
2689 "Note that activating steepest-descent energy minimization via the "
2690 "integrator .mdp option and the command gmx mdrun may "
2691 "be available in a different form in a future version of GROMACS, "
2692 "e.g. gmx minimize and an .mdp option.");
2694 /* Create 2 states on the stack and extract pointers that we will swap */
2695 em_state_t s0{}, s1{};
2696 em_state_t* s_min = &s0;
2697 em_state_t* s_try = &s1;
2699 /* Init em and store the local state in s_try */
2718 const bool simulationsShareState = false;
2719 gmx_mdoutf* outf = init_mdoutf(fplog,
2730 StartingBehavior::NewSimulation,
2731 simulationsShareState,
2733 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf),
2739 StartingBehavior::NewSimulation,
2740 simulationsShareState,
2743 /* Print to log file */
2744 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2746 /* Set variables for stepsize (in nm). This is the largest
2747 * step that we are going to make in any direction.
2749 ustep = inputrec->em_stepsize;
2752 /* Max number of steps */
2753 nsteps = inputrec->nsteps;
2757 /* Print to the screen */
2758 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2762 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2764 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
2765 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2766 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
2768 /**** HERE STARTS THE LOOP ****
2769 * count is the counter for the number of steps
2770 * bDone will be TRUE when the minimization has converged
2771 * bAbort will be TRUE when nsteps steps have been performed or when
2772 * the stepsize becomes smaller than is reasonable for machine precision
2777 while (!bDone && !bAbort)
2779 bAbort = (nsteps >= 0) && (count == nsteps);
2781 /* set new coordinates, except for first step */
2782 bool validStep = true;
2785 validStep = do_em_step(
2786 cr, inputrec, mdatoms, s_min, stepsize, s_min->f.view().forceWithPadding(), s_try, constr, count);
2791 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2795 // Signal constraint error during stepping with energy=inf
2796 s_try->epot = std::numeric_limits<real>::infinity();
2801 EnergyOutput::printHeader(fplog, count, count);
2806 s_min->epot = s_try->epot;
2809 /* Print it if necessary */
2812 if (mdrunOptions.verbose)
2815 "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2821 ((count == 0) || (s_try->epot < s_min->epot)) ? '\n' : '\r');
2825 if ((count == 0) || (s_try->epot < s_min->epot))
2827 /* Store the new (lower) energies */
2828 matrix nullBox = {};
2829 energyOutput.addDataAtEnergyStep(false,
2831 static_cast<double>(count),
2847 imdSession->fillEnergyRecord(count, TRUE);
2849 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2850 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2851 energyOutput.printStepToEnergyFile(
2852 mdoutf_get_fp_ene(outf), TRUE, do_dr, do_or, fplog, count, count, fr->fcdata.get(), nullptr);
2857 /* Now if the new energy is smaller than the previous...
2858 * or if this is the first step!
2859 * or if we did random steps!
2862 if ((count == 0) || (s_try->epot < s_min->epot))
2866 /* Test whether the convergence criterion is met... */
2867 bDone = (s_try->fmax < inputrec->em_tol);
2869 /* Copy the arrays for force, positions and energy */
2870 /* The 'Min' array always holds the coords and forces of the minimal
2872 swap_em_state(&s_min, &s_try);
2878 /* Write to trn, if necessary */
2879 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2880 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2882 fplog, cr, outf, do_x, do_f, nullptr, top_global, inputrec, count, s_min, state_global, observablesHistory);
2886 /* If energy is not smaller make the step smaller... */
2889 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2891 /* Reload the old state */
2892 em_dd_partition_system(fplog,
2911 // If the force is very small after finishing minimization,
2912 // we risk dividing by zero when calculating the step size.
2913 // So we check first if the minimization has stopped before
2914 // trying to obtain a new step size.
2917 /* Determine new step */
2918 stepsize = ustep / s_min->fmax;
2921 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2923 if (count == nsteps || ustep < 1e-12)
2925 if (count == nsteps || ustep < 1e-6)
2930 warn_step(fplog, inputrec->em_tol, s_min->fmax, count == nsteps, constr != nullptr);
2935 /* Send IMD energies and positions, if bIMD is TRUE. */
2936 if (imdSession->run(count,
2938 MASTER(cr) ? state_global->box : nullptr,
2939 MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>(),
2943 imdSession->sendPositionsAndEnergies();
2947 } /* End of the loop */
2949 /* Print some data... */
2952 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2954 write_em_traj(fplog,
2958 inputrec->nstfout != 0,
2959 ftp2fn(efSTO, nfile, fnm),
2965 observablesHistory);
2969 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2971 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2972 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps, s_min, sqrtNumAtoms);
2975 finish_em(cr, outf, walltime_accounting, wcycle);
2977 /* To print the actual number of steps we needed somewhere */
2979 // TODO: Avoid changing inputrec (#3854)
2980 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
2981 nonConstInputrec->nsteps = count;
2984 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2987 void LegacySimulator::do_nm()
2989 const char* NM = "Normal Mode Analysis";
2991 gmx_localtop_t top(top_global.ffparams);
2992 gmx_global_stat_t gstat;
2994 rvec mu_tot = { 0 };
2996 gmx_bool bSparse; /* use sparse matrix storage format */
2998 gmx_sparsematrix_t* sparse_matrix = nullptr;
2999 real* full_matrix = nullptr;
3001 /* added with respect to mdrun */
3003 real der_range = 10.0 * std::sqrt(GMX_REAL_EPS);
3005 bool bIsMaster = MASTER(cr);
3006 auto mdatoms = mdAtoms->mdatoms();
3011 "Note that activating normal-mode analysis via the integrator "
3012 ".mdp option and the command gmx mdrun may "
3013 "be available in a different form in a future version of GROMACS, "
3014 "e.g. gmx normal-modes.");
3016 if (constr != nullptr)
3020 "Constraints present with Normal Mode Analysis, this combination is not supported");
3023 gmx_shellfc_t* shellfc;
3025 em_state_t state_work{};
3027 /* Init em and store the local state in state_minimum */
3046 const bool simulationsShareState = false;
3047 gmx_mdoutf* outf = init_mdoutf(fplog,
3058 StartingBehavior::NewSimulation,
3059 simulationsShareState,
3062 std::vector<int> atom_index = get_atom_index(top_global);
3063 std::vector<gmx::RVec> fneg(atom_index.size(), { 0, 0, 0 });
3064 snew(dfdx, atom_index.size());
3070 "NOTE: This version of GROMACS has been compiled in single precision,\n"
3071 " which MIGHT not be accurate enough for normal mode analysis.\n"
3072 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
3073 " are fairly modest even if you recompile in double precision.\n\n");
3077 /* Check if we can/should use sparse storage format.
3079 * Sparse format is only useful when the Hessian itself is sparse, which it
3080 * will be when we use a cutoff.
3081 * For small systems (n<1000) it is easier to always use full matrix format, though.
3083 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
3085 GMX_LOG(mdlog.warning)
3086 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
3089 else if (atom_index.size() < 1000)
3091 GMX_LOG(mdlog.warning)
3092 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
3098 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
3102 /* Number of dimensions, based on real atoms, that is not vsites or shell */
3103 sz = DIM * atom_index.size();
3105 fprintf(stderr, "Allocating Hessian memory...\n\n");
3109 sparse_matrix = gmx_sparsematrix_init(sz);
3110 sparse_matrix->compressed_symmetric = TRUE;
3114 snew(full_matrix, sz * sz);
3117 /* Write start time and temperature */
3118 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
3120 /* fudge nr of steps to nr of atoms */
3122 // TODO: Avoid changing inputrec (#3854)
3123 auto* nonConstInputrec = const_cast<t_inputrec*>(inputrec);
3124 nonConstInputrec->nsteps = atom_index.size() * 2;
3130 "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
3135 nnodes = cr->nnodes;
3137 /* Make evaluate_energy do a single node force calculation */
3139 EnergyEvaluator energyEvaluator{ fplog, mdlog, cr, ms, top_global, &top,
3140 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
3141 vsite, constr, mdAtoms, fr, runScheduleWork, enerd };
3142 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
3143 cr->nnodes = nnodes;
3145 /* if forces are not small, warn user */
3146 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
3148 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
3149 if (state_work.fmax > 1.0e-3)
3151 GMX_LOG(mdlog.warning)
3153 "The force is probably not small enough to "
3154 "ensure that you are at a minimum.\n"
3155 "Be aware that negative eigenvalues may occur\n"
3156 "when the resulting matrix is diagonalized.");
3159 /***********************************************************
3161 * Loop over all pairs in matrix
3163 * do_force called twice. Once with positive and
3164 * once with negative displacement
3166 ************************************************************/
3168 /* Steps are divided one by one over the nodes */
3170 auto state_work_x = makeArrayRef(state_work.s.x);
3171 auto state_work_f = state_work.f.view().force();
3172 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
3174 size_t atom = atom_index[aid];
3175 for (size_t d = 0; d < DIM; d++)
3178 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
3181 x_min = state_work_x[atom][d];
3183 for (unsigned int dx = 0; (dx < 2); dx++)
3187 state_work_x[atom][d] = x_min - der_range;
3191 state_work_x[atom][d] = x_min + der_range;
3194 /* Make evaluate_energy do a single node force calculation */
3198 /* Now is the time to relax the shells */
3199 relax_shell_flexcon(fplog,
3202 mdrunOptions.verbose,
3213 state_work.s.natoms,
3214 state_work.s.x.arrayRefWithPadding(),
3215 state_work.s.v.arrayRefWithPadding(),
3217 state_work.s.lambda,
3219 &state_work.f.view(),
3230 DDBalanceRegionHandler(nullptr));
3236 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid * 2 + dx, FALSE);
3239 cr->nnodes = nnodes;
3243 std::copy(state_work_f.begin(), state_work_f.begin() + atom_index.size(), fneg.begin());
3247 /* x is restored to original */
3248 state_work_x[atom][d] = x_min;
3250 for (size_t j = 0; j < atom_index.size(); j++)
3252 for (size_t k = 0; (k < DIM); k++)
3254 dfdx[j][k] = -(state_work_f[atom_index[j]][k] - fneg[j][k]) / (2 * der_range);
3261 # define mpi_type GMX_MPI_REAL
3262 MPI_Send(dfdx[0], atom_index.size() * DIM, mpi_type, MASTER(cr), cr->nodeid, cr->mpi_comm_mygroup);
3267 for (index node = 0; (node < nnodes && aid + node < ssize(atom_index)); node++)
3273 MPI_Recv(dfdx[0], atom_index.size() * DIM, mpi_type, node, node, cr->mpi_comm_mygroup, &stat);
3278 row = (aid + node) * DIM + d;
3280 for (size_t j = 0; j < atom_index.size(); j++)
3282 for (size_t k = 0; k < DIM; k++)
3288 if (col >= row && dfdx[j][k] != 0.0)
3290 gmx_sparsematrix_increment_value(sparse_matrix, row, col, dfdx[j][k]);
3295 full_matrix[row * sz + col] = dfdx[j][k];
3302 if (mdrunOptions.verbose && fplog)
3307 /* write progress */
3308 if (bIsMaster && mdrunOptions.verbose)
3311 "\rFinished step %d out of %td",
3312 std::min<int>(atom + nnodes, atom_index.size()),
3320 fprintf(stderr, "\n\nWriting Hessian...\n");
3321 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
3324 finish_em(cr, outf, walltime_accounting, wcycle);
3326 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size() * 2);