2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file defines integrators for energy minimization
41 * \author Berk Hess <hess@kth.se>
42 * \author Erik Lindahl <erik@kth.se>
43 * \ingroup module_mdrun
56 #include "gromacs/commandline/filenm.h"
57 #include "gromacs/domdec/collect.h"
58 #include "gromacs/domdec/dlbtiming.h"
59 #include "gromacs/domdec/domdec.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/domdec/mdsetup.h"
62 #include "gromacs/domdec/partition.h"
63 #include "gromacs/ewald/pme.h"
64 #include "gromacs/fileio/confio.h"
65 #include "gromacs/fileio/mtxio.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gmxlib/nrnb.h"
68 #include "gromacs/imd/imd.h"
69 #include "gromacs/linearalgebra/sparsematrix.h"
70 #include "gromacs/listed_forces/manage_threading.h"
71 #include "gromacs/math/functions.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/mdlib/constr.h"
74 #include "gromacs/mdlib/dispersioncorrection.h"
75 #include "gromacs/mdlib/ebin.h"
76 #include "gromacs/mdlib/enerdata_utils.h"
77 #include "gromacs/mdlib/energyoutput.h"
78 #include "gromacs/mdlib/force.h"
79 #include "gromacs/mdlib/forcerec.h"
80 #include "gromacs/mdlib/gmx_omp_nthreads.h"
81 #include "gromacs/mdlib/md_support.h"
82 #include "gromacs/mdlib/mdatoms.h"
83 #include "gromacs/mdlib/stat.h"
84 #include "gromacs/mdlib/tgroup.h"
85 #include "gromacs/mdlib/trajectory_writing.h"
86 #include "gromacs/mdlib/update.h"
87 #include "gromacs/mdlib/vsite.h"
88 #include "gromacs/mdrunutility/handlerestart.h"
89 #include "gromacs/mdrunutility/printtime.h"
90 #include "gromacs/mdtypes/commrec.h"
91 #include "gromacs/mdtypes/inputrec.h"
92 #include "gromacs/mdtypes/md_enums.h"
93 #include "gromacs/mdtypes/mdrunoptions.h"
94 #include "gromacs/mdtypes/state.h"
95 #include "gromacs/pbcutil/mshift.h"
96 #include "gromacs/pbcutil/pbc.h"
97 #include "gromacs/timing/wallcycle.h"
98 #include "gromacs/timing/walltime_accounting.h"
99 #include "gromacs/topology/mtop_util.h"
100 #include "gromacs/topology/topology.h"
101 #include "gromacs/utility/cstringutil.h"
102 #include "gromacs/utility/exceptions.h"
103 #include "gromacs/utility/fatalerror.h"
104 #include "gromacs/utility/logger.h"
105 #include "gromacs/utility/smalloc.h"
107 #include "legacysimulator.h"
110 using gmx::MdrunScheduleWorkload;
112 //! Utility structure for manipulating states during EM
114 //! Copy of the global state
117 PaddedHostVector<gmx::RVec> f;
120 //! Norm of the force
128 //! Print the EM starting conditions
129 static void print_em_start(FILE *fplog,
131 gmx_walltime_accounting_t walltime_accounting,
132 gmx_wallcycle_t wcycle,
135 walltime_accounting_start_time(walltime_accounting);
136 wallcycle_start(wcycle, ewcRUN);
137 print_start(fplog, cr, walltime_accounting, name);
140 //! Stop counting time for EM
141 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
142 gmx_wallcycle_t wcycle)
144 wallcycle_stop(wcycle, ewcRUN);
146 walltime_accounting_end_time(walltime_accounting);
149 //! Printing a log file and console header
150 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
153 fprintf(out, "%s:\n", minimizer);
154 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
155 fprintf(out, " Number of steps = %12d\n", nsteps);
158 //! Print warning message
159 static void warn_step(FILE *fp,
165 constexpr bool realIsDouble = GMX_DOUBLE;
168 if (!std::isfinite(fmax))
171 "\nEnergy minimization has stopped because the force "
172 "on at least one atom is not finite. This usually means "
173 "atoms are overlapping. Modify the input coordinates to "
174 "remove atom overlap or use soft-core potentials with "
175 "the free energy code to avoid infinite forces.\n%s",
177 "You could also be lucky that switching to double precision "
178 "is sufficient to obtain finite forces.\n" :
184 "\nEnergy minimization reached the maximum number "
185 "of steps before the forces reached the requested "
186 "precision Fmax < %g.\n", ftol);
191 "\nEnergy minimization has stopped, but the forces have "
192 "not converged to the requested precision Fmax < %g (which "
193 "may not be possible for your system). It stopped "
194 "because the algorithm tried to make a new step whose size "
195 "was too small, or there was no change in the energy since "
196 "last step. Either way, we regard the minimization as "
197 "converged to within the available machine precision, "
198 "given your starting configuration and EM parameters.\n%s%s",
201 "\nDouble precision normally gives you higher accuracy, but "
202 "this is often not needed for preparing to run molecular "
206 "You might need to increase your constraint accuracy, or turn\n"
207 "off constraints altogether (set constraints = none in mdp file)\n" :
211 fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
212 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
215 //! Print message about convergence of the EM
216 static void print_converged(FILE *fp, const char *alg, real ftol,
217 int64_t count, gmx_bool bDone, int64_t nsteps,
218 const em_state_t *ems, double sqrtNumAtoms)
220 char buf[STEPSTRSIZE];
224 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
225 alg, ftol, gmx_step_str(count, buf));
227 else if (count < nsteps)
229 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
230 "but did not reach the requested Fmax < %g.\n",
231 alg, gmx_step_str(count, buf), ftol);
235 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
236 alg, ftol, gmx_step_str(count, buf));
240 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
241 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
242 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm/sqrtNumAtoms);
244 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
245 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
246 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm/sqrtNumAtoms);
250 //! Compute the norm and max of the force array in parallel
251 static void get_f_norm_max(const t_commrec *cr,
252 t_grpopts *opts, t_mdatoms *mdatoms, const rvec *f,
253 real *fnorm, real *fmax, int *a_fmax)
257 int la_max, a_max, start, end, i, m, gf;
259 /* This routine finds the largest force and returns it.
260 * On parallel machines the global max is taken.
266 end = mdatoms->homenr;
267 if (mdatoms->cFREEZE)
269 for (i = start; i < end; i++)
271 gf = mdatoms->cFREEZE[i];
273 for (m = 0; m < DIM; m++)
275 if (!opts->nFreeze[gf][m])
277 fam += gmx::square(f[i][m]);
290 for (i = start; i < end; i++)
302 if (la_max >= 0 && DOMAINDECOMP(cr))
304 a_max = cr->dd->globalAtomIndices[la_max];
312 snew(sum, 2*cr->nnodes+1);
313 sum[2*cr->nodeid] = fmax2;
314 sum[2*cr->nodeid+1] = a_max;
315 sum[2*cr->nnodes] = fnorm2;
316 gmx_sumd(2*cr->nnodes+1, sum, cr);
317 fnorm2 = sum[2*cr->nnodes];
318 /* Determine the global maximum */
319 for (i = 0; i < cr->nnodes; i++)
321 if (sum[2*i] > fmax2)
324 a_max = gmx::roundToInt(sum[2*i+1]);
332 *fnorm = sqrt(fnorm2);
344 //! Compute the norm of the force
345 static void get_state_f_norm_max(const t_commrec *cr,
346 t_grpopts *opts, t_mdatoms *mdatoms,
349 get_f_norm_max(cr, opts, mdatoms, ems->f.rvec_array(),
350 &ems->fnorm, &ems->fmax, &ems->a_fmax);
353 //! Initialize the energy minimization
354 static void init_em(FILE *fplog,
355 const gmx::MDLogger &mdlog,
359 gmx::ImdSession *imdSession,
361 t_state *state_global, gmx_mtop_t *top_global,
362 em_state_t *ems, gmx_localtop_t *top,
365 t_graph **graph, gmx::MDAtoms *mdAtoms, gmx_global_stat_t *gstat,
366 gmx_vsite_t *vsite, gmx::Constraints *constr, gmx_shellfc_t **shellfc)
372 fprintf(fplog, "Initiating %s\n", title);
377 state_global->ngtc = 0;
379 initialize_lambdas(fplog, *ir, MASTER(cr), &(state_global->fep_state), state_global->lambda, nullptr);
383 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
385 *shellfc = init_shell_flexcon(stdout,
387 constr ? constr->numFlexibleConstraints() : 0,
393 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
395 /* With energy minimization, shells and flexible constraints are
396 * automatically minimized when treated like normal DOFS.
398 if (shellfc != nullptr)
404 auto mdatoms = mdAtoms->mdatoms();
405 if (DOMAINDECOMP(cr))
407 top->useInDomainDecomp_ = true;
408 dd_init_local_top(*top_global, top);
410 dd_init_local_state(cr->dd, state_global, &ems->s);
412 /* Distribute the charge groups over the nodes from the master node */
413 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
414 state_global, *top_global, ir, imdSession, pull_work,
415 &ems->s, &ems->f, mdAtoms, top,
417 nrnb, nullptr, FALSE);
418 dd_store_state(cr->dd, &ems->s);
424 state_change_natoms(state_global, state_global->natoms);
425 /* Just copy the state */
426 ems->s = *state_global;
427 state_change_natoms(&ems->s, ems->s.natoms);
428 ems->f.resizeWithPadding(ems->s.natoms);
430 mdAlgorithmsSetupAtomData(cr, ir, *top_global, top, fr,
432 constr, vsite, shellfc ? *shellfc : nullptr);
436 set_vsite_top(vsite, top, mdatoms);
440 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
444 // TODO how should this cross-module support dependency be managed?
445 if (ir->eConstrAlg == econtSHAKE &&
446 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
448 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
449 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
452 if (!ir->bContinuation)
454 /* Constrain the starting coordinates */
456 constr->apply(TRUE, TRUE,
458 ems->s.x.rvec_array(),
459 ems->s.x.rvec_array(),
462 ems->s.lambda[efptFEP], &dvdl_constr,
463 nullptr, nullptr, gmx::ConstraintVariable::Positions);
469 *gstat = global_stat_init(ir);
476 calc_shifts(ems->s.box, fr->shift_vec);
479 //! Finalize the minimization
480 static void finish_em(const t_commrec *cr, gmx_mdoutf_t outf,
481 gmx_walltime_accounting_t walltime_accounting,
482 gmx_wallcycle_t wcycle)
484 if (!thisRankHasDuty(cr, DUTY_PME))
486 /* Tell the PME only node to finish */
487 gmx_pme_send_finish(cr);
492 em_time_end(walltime_accounting, wcycle);
495 //! Swap two different EM states during minimization
496 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
505 //! Save the EM trajectory
506 static void write_em_traj(FILE *fplog, const t_commrec *cr,
508 gmx_bool bX, gmx_bool bF, const char *confout,
509 gmx_mtop_t *top_global,
510 t_inputrec *ir, int64_t step,
512 t_state *state_global,
513 ObservablesHistory *observablesHistory)
519 mdof_flags |= MDOF_X;
523 mdof_flags |= MDOF_F;
526 /* If we want IMD output, set appropriate MDOF flag */
529 mdof_flags |= MDOF_IMD;
532 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
533 top_global->natoms, step, static_cast<double>(step),
534 &state->s, state_global, observablesHistory,
537 if (confout != nullptr)
539 if (DOMAINDECOMP(cr))
541 /* If bX=true, x was collected to state_global in the call above */
544 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
545 dd_collect_vec(cr->dd, &state->s, state->s.x, globalXRef);
550 /* Copy the local state pointer */
551 state_global = &state->s;
556 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
558 /* Make molecules whole only for confout writing */
559 do_pbc_mtop(ir->ePBC, state->s.box, top_global,
560 state_global->x.rvec_array());
563 write_sto_conf_mtop(confout,
564 *top_global->name, top_global,
565 state_global->x.rvec_array(), nullptr, ir->ePBC, state->s.box);
570 //! \brief Do one minimization step
572 // \returns true when the step succeeded, false when a constraint error occurred
573 static bool do_em_step(const t_commrec *cr,
574 t_inputrec *ir, t_mdatoms *md,
575 em_state_t *ems1, real a, const PaddedHostVector<gmx::RVec> *force,
577 gmx::Constraints *constr,
584 int nthreads gmx_unused;
586 bool validStep = true;
591 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
593 gmx_incons("state mismatch in do_em_step");
596 s2->flags = s1->flags;
598 if (s2->natoms != s1->natoms)
600 state_change_natoms(s2, s1->natoms);
601 ems2->f.resizeWithPadding(s2->natoms);
603 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
605 s2->cg_gl.resize(s1->cg_gl.size());
608 copy_mat(s1->box, s2->box);
609 /* Copy free energy state */
610 s2->lambda = s1->lambda;
611 copy_mat(s1->box, s2->box);
616 nthreads = gmx_omp_nthreads_get(emntUpdate);
617 #pragma omp parallel num_threads(nthreads)
619 const rvec *x1 = s1->x.rvec_array();
620 rvec *x2 = s2->x.rvec_array();
621 const rvec *f = force->rvec_array();
624 #pragma omp for schedule(static) nowait
625 for (int i = start; i < end; i++)
633 for (int m = 0; m < DIM; m++)
635 if (ir->opts.nFreeze[gf][m])
641 x2[i][m] = x1[i][m] + a*f[i][m];
645 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
648 if (s2->flags & (1<<estCGP))
650 /* Copy the CG p vector */
651 const rvec *p1 = s1->cg_p.rvec_array();
652 rvec *p2 = s2->cg_p.rvec_array();
653 #pragma omp for schedule(static) nowait
654 for (int i = start; i < end; i++)
656 // Trivial OpenMP block that does not throw
657 copy_rvec(p1[i], p2[i]);
661 if (DOMAINDECOMP(cr))
663 s2->ddp_count = s1->ddp_count;
665 /* OpenMP does not supported unsigned loop variables */
666 #pragma omp for schedule(static) nowait
667 for (gmx::index i = 0; i < gmx::ssize(s2->cg_gl); i++)
669 s2->cg_gl[i] = s1->cg_gl[i];
671 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
679 constr->apply(TRUE, TRUE,
681 s1->x.rvec_array(), s2->x.rvec_array(),
683 s2->lambda[efptBONDED], &dvdl_constr,
684 nullptr, nullptr, gmx::ConstraintVariable::Positions);
688 /* This global reduction will affect performance at high
689 * parallelization, but we can not really avoid it.
690 * But usually EM is not run at high parallelization.
692 int reductionBuffer = static_cast<int>(!validStep);
693 gmx_sumi(1, &reductionBuffer, cr);
694 validStep = (reductionBuffer == 0);
697 // We should move this check to the different minimizers
698 if (!validStep && ir->eI != eiSteep)
700 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
701 EI(ir->eI), EI(eiSteep), EI(ir->eI));
708 //! Prepare EM for using domain decomposition parallellization
709 static void em_dd_partition_system(FILE *fplog,
710 const gmx::MDLogger &mdlog,
711 int step, const t_commrec *cr,
712 gmx_mtop_t *top_global, t_inputrec *ir,
713 gmx::ImdSession *imdSession,
715 em_state_t *ems, gmx_localtop_t *top,
716 gmx::MDAtoms *mdAtoms, t_forcerec *fr,
717 gmx_vsite_t *vsite, gmx::Constraints *constr,
718 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
720 /* Repartition the domain decomposition */
721 dd_partition_system(fplog, mdlog, step, cr, FALSE, 1,
722 nullptr, *top_global, ir, imdSession, pull_work,
724 mdAtoms, top, fr, vsite, constr,
725 nrnb, wcycle, FALSE);
726 dd_store_state(cr->dd, &ems->s);
732 /*! \brief Class to handle the work of setting and doing an energy evaluation.
734 * This class is a mere aggregate of parameters to pass to evaluate an
735 * energy, so that future changes to names and types of them consume
736 * less time when refactoring other code.
738 * Aggregate initialization is used, for which the chief risk is that
739 * if a member is added at the end and not all initializer lists are
740 * updated, then the member will be value initialized, which will
741 * typically mean initialization to zero.
743 * Use a braced initializer list to construct one of these. */
744 class EnergyEvaluator
747 /*! \brief Evaluates an energy on the state in \c ems.
749 * \todo In practice, the same objects mu_tot, vir, and pres
750 * are always passed to this function, so we would rather have
751 * them as data members. However, their C-array types are
752 * unsuited for aggregate initialization. When the types
753 * improve, the call signature of this method can be reduced.
755 void run(em_state_t *ems, rvec mu_tot,
756 tensor vir, tensor pres,
757 int64_t count, gmx_bool bFirst);
758 //! Handles logging (deprecated).
761 const gmx::MDLogger &mdlog;
762 //! Handles communication.
764 //! Coordinates multi-simulations.
765 const gmx_multisim_t *ms;
766 //! Holds the simulation topology.
767 gmx_mtop_t *top_global;
768 //! Holds the domain topology.
770 //! User input options.
771 t_inputrec *inputrec;
772 //! The Interactive Molecular Dynamics session.
773 gmx::ImdSession *imdSession;
774 //! The pull work object.
776 //! Manages flop accounting.
778 //! Manages wall cycle accounting.
779 gmx_wallcycle_t wcycle;
780 //! Coordinates global reduction.
781 gmx_global_stat_t gstat;
782 //! Handles virtual sites.
784 //! Handles constraints.
785 gmx::Constraints *constr;
786 //! Handles strange things.
788 //! Molecular graph for SHAKE.
790 //! Per-atom data for this domain.
791 gmx::MDAtoms *mdAtoms;
792 //! Handles how to calculate the forces.
794 //! Schedule of force-calculation work each step for this task.
795 MdrunScheduleWorkload *runScheduleWork;
796 //! Stores the computed energies.
797 gmx_enerdata_t *enerd;
801 EnergyEvaluator::run(em_state_t *ems, rvec mu_tot,
802 tensor vir, tensor pres,
803 int64_t count, gmx_bool bFirst)
807 tensor force_vir, shake_vir, ekin;
811 /* Set the time to the initial time, the time does not change during EM */
812 t = inputrec->init_t;
815 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
817 /* This is the first state or an old state used before the last ns */
823 if (inputrec->nstlist > 0)
831 construct_vsites(vsite, ems->s.x.rvec_array(), 1, nullptr,
832 top->idef.iparams, top->idef.il,
833 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
836 if (DOMAINDECOMP(cr) && bNS)
838 /* Repartition the domain decomposition */
839 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession,
841 ems, top, mdAtoms, fr, vsite, constr,
845 /* Calc force & energy on new trial position */
846 /* do_force always puts the charge groups in the box and shifts again
847 * We do not unshift, so molecules are always whole in congrad.c
849 do_force(fplog, cr, ms, inputrec, nullptr, nullptr, imdSession,
851 count, nrnb, wcycle, top,
852 ems->s.box, ems->s.x.arrayRefWithPadding(), &ems->s.hist,
853 ems->f.arrayRefWithPadding(), force_vir, mdAtoms->mdatoms(), enerd, fcd,
854 ems->s.lambda, graph, fr, runScheduleWork, vsite, mu_tot, t, nullptr,
855 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
856 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
857 (bNS ? GMX_FORCE_NS : 0),
858 DDBalanceRegionHandler(cr));
860 /* Clear the unused shake virial and pressure */
861 clear_mat(shake_vir);
864 /* Communicate stuff when parallel */
865 if (PAR(cr) && inputrec->eI != eiNM)
867 wallcycle_start(wcycle, ewcMoveE);
869 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
870 inputrec, nullptr, nullptr, nullptr, 1, &terminate,
876 wallcycle_stop(wcycle, ewcMoveE);
879 if (fr->dispersionCorrection)
881 /* Calculate long range corrections to pressure and energy */
882 const DispersionCorrection::Correction correction =
883 fr->dispersionCorrection->calculate(ems->s.box, ems->s.lambda[efptVDW]);
885 enerd->term[F_DISPCORR] = correction.energy;
886 enerd->term[F_EPOT] += correction.energy;
887 enerd->term[F_PRES] += correction.pressure;
888 enerd->term[F_DVDL] += correction.dvdl;
892 enerd->term[F_DISPCORR] = 0;
895 ems->epot = enerd->term[F_EPOT];
899 /* Project out the constraint components of the force */
901 rvec *f_rvec = ems->f.rvec_array();
902 constr->apply(FALSE, FALSE,
904 ems->s.x.rvec_array(), f_rvec, f_rvec,
906 ems->s.lambda[efptBONDED], &dvdl_constr,
907 nullptr, &shake_vir, gmx::ConstraintVariable::ForceDispl);
908 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
909 m_add(force_vir, shake_vir, vir);
913 copy_mat(force_vir, vir);
917 enerd->term[F_PRES] =
918 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
920 sum_dhdl(enerd, ems->s.lambda, *inputrec->fepvals);
922 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
924 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
930 //! Parallel utility summing energies and forces
931 static double reorder_partsum(const t_commrec *cr, t_grpopts *opts,
932 gmx_mtop_t *top_global,
933 em_state_t *s_min, em_state_t *s_b)
937 fprintf(debug, "Doing reorder_partsum\n");
940 const rvec *fm = s_min->f.rvec_array();
941 const rvec *fb = s_b->f.rvec_array();
943 /* Collect fm in a global vector fmg.
944 * This conflicts with the spirit of domain decomposition,
945 * but to fully optimize this a much more complicated algorithm is required.
947 const int natoms = top_global->natoms;
951 gmx::ArrayRef<const int> indicesMin = s_b->s.cg_gl;
953 for (int a : indicesMin)
955 copy_rvec(fm[i], fmg[a]);
957 gmx_sum(top_global->natoms*3, fmg[0], cr);
959 /* Now we will determine the part of the sum for the cgs in state s_b */
960 gmx::ArrayRef<const int> indicesB = s_b->s.cg_gl;
965 gmx::ArrayRef<unsigned char> grpnrFREEZE = top_global->groups.groupNumbers[SimulationAtomGroupType::Freeze];
966 for (int a : indicesB)
968 if (!grpnrFREEZE.empty())
972 for (int m = 0; m < DIM; m++)
974 if (!opts->nFreeze[gf][m])
976 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
987 //! Print some stuff, like beta, whatever that means.
988 static real pr_beta(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
989 gmx_mtop_t *top_global,
990 em_state_t *s_min, em_state_t *s_b)
994 /* This is just the classical Polak-Ribiere calculation of beta;
995 * it looks a bit complicated since we take freeze groups into account,
996 * and might have to sum it in parallel runs.
999 if (!DOMAINDECOMP(cr) ||
1000 (s_min->s.ddp_count == cr->dd->ddp_count &&
1001 s_b->s.ddp_count == cr->dd->ddp_count))
1003 const rvec *fm = s_min->f.rvec_array();
1004 const rvec *fb = s_b->f.rvec_array();
1007 /* This part of code can be incorrect with DD,
1008 * since the atom ordering in s_b and s_min might differ.
1010 for (int i = 0; i < mdatoms->homenr; i++)
1012 if (mdatoms->cFREEZE)
1014 gf = mdatoms->cFREEZE[i];
1016 for (int m = 0; m < DIM; m++)
1018 if (!opts->nFreeze[gf][m])
1020 sum += (fb[i][m] - fm[i][m])*fb[i][m];
1027 /* We need to reorder cgs while summing */
1028 sum = reorder_partsum(cr, opts, top_global, s_min, s_b);
1032 gmx_sumd(1, &sum, cr);
1035 return sum/gmx::square(s_min->fnorm);
1042 LegacySimulator::do_cg()
1044 const char *CG = "Polak-Ribiere Conjugate Gradients";
1047 gmx_global_stat_t gstat;
1049 double tmp, minstep;
1051 real a, b, c, beta = 0.0;
1054 gmx_bool converged, foundlower;
1056 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1058 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1059 int m, step, nminstep;
1060 auto mdatoms = mdAtoms->mdatoms();
1062 GMX_LOG(mdlog.info).asParagraph().
1063 appendText("Note that activating conjugate gradient energy minimization via the "
1064 "integrator .mdp option and the command gmx mdrun may "
1065 "be available in a different form in a future version of GROMACS, "
1066 "e.g. gmx minimize and an .mdp option.");
1072 // In CG, the state is extended with a search direction
1073 state_global->flags |= (1<<estCGP);
1075 // Ensure the extra per-atom state array gets allocated
1076 state_change_natoms(state_global, state_global->natoms);
1078 // Initialize the search direction to zero
1079 for (RVec &cg_p : state_global->cg_p)
1085 /* Create 4 states on the stack and extract pointers that we will swap */
1086 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1087 em_state_t *s_min = &s0;
1088 em_state_t *s_a = &s1;
1089 em_state_t *s_b = &s2;
1090 em_state_t *s_c = &s3;
1092 /* Init em and store the local state in s_min */
1093 init_em(fplog, mdlog, CG, cr, inputrec, imdSession,
1095 state_global, top_global, s_min, &top,
1096 nrnb, fr, &graph, mdAtoms, &gstat,
1097 vsite, constr, nullptr);
1098 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1099 StartingBehavior::NewSimulation);
1100 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr, false, mdModulesNotifier);
1102 /* Print to log file */
1103 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1105 /* Max number of steps */
1106 number_steps = inputrec->nsteps;
1110 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1114 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1117 EnergyEvaluator energyEvaluator {
1118 fplog, mdlog, cr, ms,
1120 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1121 vsite, constr, fcd, graph,
1122 mdAtoms, fr, runScheduleWork, enerd
1124 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1125 /* do_force always puts the charge groups in the box and shifts again
1126 * We do not unshift, so molecules are always whole in congrad.c
1128 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1132 /* Copy stuff to the energy bin for easy printing etc. */
1133 matrix nullBox = {};
1134 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1135 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1136 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1138 energyOutput.printHeader(fplog, step, step);
1139 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE,
1144 /* Estimate/guess the initial stepsize */
1145 stepsize = inputrec->em_stepsize/s_min->fnorm;
1149 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1150 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1151 s_min->fmax, s_min->a_fmax+1);
1152 fprintf(stderr, " F-Norm = %12.5e\n",
1153 s_min->fnorm/sqrtNumAtoms);
1154 fprintf(stderr, "\n");
1155 /* and copy to the log file too... */
1156 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1157 s_min->fmax, s_min->a_fmax+1);
1158 fprintf(fplog, " F-Norm = %12.5e\n",
1159 s_min->fnorm/sqrtNumAtoms);
1160 fprintf(fplog, "\n");
1162 /* Start the loop over CG steps.
1163 * Each successful step is counted, and we continue until
1164 * we either converge or reach the max number of steps.
1167 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1170 /* start taking steps in a new direction
1171 * First time we enter the routine, beta=0, and the direction is
1172 * simply the negative gradient.
1175 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1176 rvec *pm = s_min->s.cg_p.rvec_array();
1177 const rvec *sfm = s_min->f.rvec_array();
1180 for (int i = 0; i < mdatoms->homenr; i++)
1182 if (mdatoms->cFREEZE)
1184 gf = mdatoms->cFREEZE[i];
1186 for (m = 0; m < DIM; m++)
1188 if (!inputrec->opts.nFreeze[gf][m])
1190 pm[i][m] = sfm[i][m] + beta*pm[i][m];
1191 gpa -= pm[i][m]*sfm[i][m];
1192 /* f is negative gradient, thus the sign */
1201 /* Sum the gradient along the line across CPUs */
1204 gmx_sumd(1, &gpa, cr);
1207 /* Calculate the norm of the search vector */
1208 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1210 /* Just in case stepsize reaches zero due to numerical precision... */
1213 stepsize = inputrec->em_stepsize/pnorm;
1217 * Double check the value of the derivative in the search direction.
1218 * If it is positive it must be due to the old information in the
1219 * CG formula, so just remove that and start over with beta=0.
1220 * This corresponds to a steepest descent step.
1225 step--; /* Don't count this step since we are restarting */
1226 continue; /* Go back to the beginning of the big for-loop */
1229 /* Calculate minimum allowed stepsize, before the average (norm)
1230 * relative change in coordinate is smaller than precision
1233 auto s_min_x = makeArrayRef(s_min->s.x);
1234 for (int i = 0; i < mdatoms->homenr; i++)
1236 for (m = 0; m < DIM; m++)
1238 tmp = fabs(s_min_x[i][m]);
1247 /* Add up from all CPUs */
1250 gmx_sumd(1, &minstep, cr);
1253 minstep = GMX_REAL_EPS/sqrt(minstep/(3*top_global->natoms));
1255 if (stepsize < minstep)
1261 /* Write coordinates if necessary */
1262 do_x = do_per_step(step, inputrec->nstxout);
1263 do_f = do_per_step(step, inputrec->nstfout);
1265 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1266 top_global, inputrec, step,
1267 s_min, state_global, observablesHistory);
1269 /* Take a step downhill.
1270 * In theory, we should minimize the function along this direction.
1271 * That is quite possible, but it turns out to take 5-10 function evaluations
1272 * for each line. However, we dont really need to find the exact minimum -
1273 * it is much better to start a new CG step in a modified direction as soon
1274 * as we are close to it. This will save a lot of energy evaluations.
1276 * In practice, we just try to take a single step.
1277 * If it worked (i.e. lowered the energy), we increase the stepsize but
1278 * the continue straight to the next CG step without trying to find any minimum.
1279 * If it didn't work (higher energy), there must be a minimum somewhere between
1280 * the old position and the new one.
1282 * Due to the finite numerical accuracy, it turns out that it is a good idea
1283 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1284 * This leads to lower final energies in the tests I've done. / Erik
1286 s_a->epot = s_min->epot;
1288 c = a + stepsize; /* reference position along line is zero */
1290 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1292 em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec, imdSession,
1294 s_min, &top, mdAtoms, fr, vsite, constr,
1298 /* Take a trial step (new coords in s_c) */
1299 do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c,
1303 /* Calculate energy for the trial step */
1304 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1306 /* Calc derivative along line */
1307 const rvec *pc = s_c->s.cg_p.rvec_array();
1308 const rvec *sfc = s_c->f.rvec_array();
1310 for (int i = 0; i < mdatoms->homenr; i++)
1312 for (m = 0; m < DIM; m++)
1314 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1317 /* Sum the gradient along the line across CPUs */
1320 gmx_sumd(1, &gpc, cr);
1323 /* This is the max amount of increase in energy we tolerate */
1324 tmp = std::sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1326 /* Accept the step if the energy is lower, or if it is not significantly higher
1327 * and the line derivative is still negative.
1329 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1332 /* Great, we found a better energy. Increase step for next iteration
1333 * if we are still going down, decrease it otherwise
1337 stepsize *= 1.618034; /* The golden section */
1341 stepsize *= 0.618034; /* 1/golden section */
1346 /* New energy is the same or higher. We will have to do some work
1347 * to find a smaller value in the interval. Take smaller step next time!
1350 stepsize *= 0.618034;
1356 /* OK, if we didn't find a lower value we will have to locate one now - there must
1357 * be one in the interval [a=0,c].
1358 * The same thing is valid here, though: Don't spend dozens of iterations to find
1359 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1360 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1362 * I also have a safeguard for potentially really pathological functions so we never
1363 * take more than 20 steps before we give up ...
1365 * If we already found a lower value we just skip this step and continue to the update.
1374 /* Select a new trial point.
1375 * If the derivatives at points a & c have different sign we interpolate to zero,
1376 * otherwise just do a bisection.
1378 if (gpa < 0 && gpc > 0)
1380 b = a + gpa*(a-c)/(gpc-gpa);
1387 /* safeguard if interpolation close to machine accuracy causes errors:
1388 * never go outside the interval
1390 if (b <= a || b >= c)
1395 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1397 /* Reload the old state */
1398 em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec, imdSession,
1400 s_min, &top, mdAtoms, fr, vsite, constr,
1404 /* Take a trial step to this new point - new coords in s_b */
1405 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b,
1409 /* Calculate energy for the trial step */
1410 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1412 /* p does not change within a step, but since the domain decomposition
1413 * might change, we have to use cg_p of s_b here.
1415 const rvec *pb = s_b->s.cg_p.rvec_array();
1416 const rvec *sfb = s_b->f.rvec_array();
1418 for (int i = 0; i < mdatoms->homenr; i++)
1420 for (m = 0; m < DIM; m++)
1422 gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1425 /* Sum the gradient along the line across CPUs */
1428 gmx_sumd(1, &gpb, cr);
1433 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1434 s_a->epot, s_b->epot, s_c->epot, gpb);
1437 epot_repl = s_b->epot;
1439 /* Keep one of the intervals based on the value of the derivative at the new point */
1442 /* Replace c endpoint with b */
1443 swap_em_state(&s_b, &s_c);
1449 /* Replace a endpoint with b */
1450 swap_em_state(&s_b, &s_a);
1456 * Stop search as soon as we find a value smaller than the endpoints.
1457 * Never run more than 20 steps, no matter what.
1461 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1464 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1467 /* OK. We couldn't find a significantly lower energy.
1468 * If beta==0 this was steepest descent, and then we give up.
1469 * If not, set beta=0 and restart with steepest descent before quitting.
1479 /* Reset memory before giving up */
1485 /* Select min energy state of A & C, put the best in B.
1487 if (s_c->epot < s_a->epot)
1491 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1492 s_c->epot, s_a->epot);
1494 swap_em_state(&s_b, &s_c);
1501 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1502 s_a->epot, s_c->epot);
1504 swap_em_state(&s_b, &s_a);
1513 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1516 swap_em_state(&s_b, &s_c);
1520 /* new search direction */
1521 /* beta = 0 means forget all memory and restart with steepest descents. */
1522 if (nstcg && ((step % nstcg) == 0))
1528 /* s_min->fnorm cannot be zero, because then we would have converged
1532 /* Polak-Ribiere update.
1533 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1535 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1537 /* Limit beta to prevent oscillations */
1538 if (fabs(beta) > 5.0)
1544 /* update positions */
1545 swap_em_state(&s_min, &s_b);
1548 /* Print it if necessary */
1551 if (mdrunOptions.verbose)
1553 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1554 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1555 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1556 s_min->fmax, s_min->a_fmax+1);
1559 /* Store the new (lower) energies */
1560 matrix nullBox = {};
1561 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1562 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1563 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1565 do_log = do_per_step(step, inputrec->nstlog);
1566 do_ene = do_per_step(step, inputrec->nstenergy);
1568 imdSession->fillEnergyRecord(step, TRUE);
1572 energyOutput.printHeader(fplog, step, step);
1574 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1575 do_log ? fplog : nullptr, step, step,
1579 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1580 if (MASTER(cr) && imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0))
1582 imdSession->sendPositionsAndEnergies();
1585 /* Stop when the maximum force lies below tolerance.
1586 * If we have reached machine precision, converged is already set to true.
1588 converged = converged || (s_min->fmax < inputrec->em_tol);
1590 } /* End of the loop */
1594 step--; /* we never took that last step in this case */
1597 if (s_min->fmax > inputrec->em_tol)
1601 warn_step(fplog, inputrec->em_tol, s_min->fmax,
1602 step-1 == number_steps, FALSE);
1609 /* If we printed energy and/or logfile last step (which was the last step)
1610 * we don't have to do it again, but otherwise print the final values.
1614 /* Write final value to log since we didn't do anything the last step */
1615 energyOutput.printHeader(fplog, step, step);
1617 if (!do_ene || !do_log)
1619 /* Write final energy file entries */
1620 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1621 !do_log ? fplog : nullptr, step, step,
1626 /* Print some stuff... */
1629 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1633 * For accurate normal mode calculation it is imperative that we
1634 * store the last conformation into the full precision binary trajectory.
1636 * However, we should only do it if we did NOT already write this step
1637 * above (which we did if do_x or do_f was true).
1639 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1640 * in the trajectory with the same step number.
1642 do_x = !do_per_step(step, inputrec->nstxout);
1643 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1645 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1646 top_global, inputrec, step,
1647 s_min, state_global, observablesHistory);
1652 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1653 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1654 s_min, sqrtNumAtoms);
1655 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1656 s_min, sqrtNumAtoms);
1658 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1661 finish_em(cr, outf, walltime_accounting, wcycle);
1663 /* To print the actual number of steps we needed somewhere */
1664 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1669 LegacySimulator::do_lbfgs()
1671 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1674 gmx_global_stat_t gstat;
1676 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1677 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1678 real *rho, *alpha, *p, *s, **dx, **dg;
1679 real a, b, c, maxdelta, delta;
1681 real dgdx, dgdg, sq, yr, beta;
1684 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1686 int start, end, number_steps;
1687 int i, k, m, n, gf, step;
1689 auto mdatoms = mdAtoms->mdatoms();
1691 GMX_LOG(mdlog.info).asParagraph().
1692 appendText("Note that activating L-BFGS energy minimization via the "
1693 "integrator .mdp option and the command gmx mdrun may "
1694 "be available in a different form in a future version of GROMACS, "
1695 "e.g. gmx minimize and an .mdp option.");
1699 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1702 if (nullptr != constr)
1704 gmx_fatal(FARGS, "The combination of constraints and L-BFGS minimization is not implemented. Either do not use constraints, or use another minimizer (e.g. steepest descent).");
1707 n = 3*state_global->natoms;
1708 nmaxcorr = inputrec->nbfgscorr;
1713 snew(rho, nmaxcorr);
1714 snew(alpha, nmaxcorr);
1717 for (i = 0; i < nmaxcorr; i++)
1723 for (i = 0; i < nmaxcorr; i++)
1732 init_em(fplog, mdlog, LBFGS, cr, inputrec, imdSession,
1734 state_global, top_global, &ems, &top,
1735 nrnb, fr, &graph, mdAtoms, &gstat,
1736 vsite, constr, nullptr);
1737 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
1738 StartingBehavior::NewSimulation);
1739 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr, false, mdModulesNotifier);
1742 end = mdatoms->homenr;
1744 /* We need 4 working states */
1745 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1746 em_state_t *sa = &s0;
1747 em_state_t *sb = &s1;
1748 em_state_t *sc = &s2;
1749 em_state_t *last = &s3;
1750 /* Initialize by copying the state from ems (we could skip x and f here) */
1755 /* Print to log file */
1756 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1758 do_log = do_ene = do_x = do_f = TRUE;
1760 /* Max number of steps */
1761 number_steps = inputrec->nsteps;
1763 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1765 for (i = start; i < end; i++)
1767 if (mdatoms->cFREEZE)
1769 gf = mdatoms->cFREEZE[i];
1771 for (m = 0; m < DIM; m++)
1773 frozen[3*i+m] = (inputrec->opts.nFreeze[gf][m] != 0);
1778 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1782 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1787 construct_vsites(vsite, state_global->x.rvec_array(), 1, nullptr,
1788 top.idef.iparams, top.idef.il,
1789 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1792 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1793 /* do_force always puts the charge groups in the box and shifts again
1794 * We do not unshift, so molecules are always whole
1797 EnergyEvaluator energyEvaluator {
1798 fplog, mdlog, cr, ms,
1800 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
1801 vsite, constr, fcd, graph,
1802 mdAtoms, fr, runScheduleWork, enerd
1804 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1808 /* Copy stuff to the energy bin for easy printing etc. */
1809 matrix nullBox = {};
1810 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
1811 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1812 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1814 energyOutput.printHeader(fplog, step, step);
1815 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE,
1820 /* Set the initial step.
1821 * since it will be multiplied by the non-normalized search direction
1822 * vector (force vector the first time), we scale it by the
1823 * norm of the force.
1828 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1829 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1830 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1831 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1832 fprintf(stderr, "\n");
1833 /* and copy to the log file too... */
1834 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1835 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1836 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1837 fprintf(fplog, "\n");
1840 // Point is an index to the memory of search directions, where 0 is the first one.
1843 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1844 real *fInit = static_cast<real *>(ems.f.rvec_array()[0]);
1845 for (i = 0; i < n; i++)
1849 dx[point][i] = fInit[i]; /* Initial search direction */
1857 // Stepsize will be modified during the search, and actually it is not critical
1858 // (the main efficiency in the algorithm comes from changing directions), but
1859 // we still need an initial value, so estimate it as the inverse of the norm
1860 // so we take small steps where the potential fluctuates a lot.
1861 stepsize = 1.0/ems.fnorm;
1863 /* Start the loop over BFGS steps.
1864 * Each successful step is counted, and we continue until
1865 * we either converge or reach the max number of steps.
1870 /* Set the gradient from the force */
1872 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1875 /* Write coordinates if necessary */
1876 do_x = do_per_step(step, inputrec->nstxout);
1877 do_f = do_per_step(step, inputrec->nstfout);
1882 mdof_flags |= MDOF_X;
1887 mdof_flags |= MDOF_F;
1892 mdof_flags |= MDOF_IMD;
1895 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1896 top_global->natoms, step, static_cast<real>(step), &ems.s,
1897 state_global, observablesHistory, ems.f);
1899 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1901 /* make s a pointer to current search direction - point=0 first time we get here */
1904 real *xx = static_cast<real *>(ems.s.x.rvec_array()[0]);
1905 real *ff = static_cast<real *>(ems.f.rvec_array()[0]);
1907 // calculate line gradient in position A
1908 for (gpa = 0, i = 0; i < n; i++)
1913 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1914 * relative change in coordinate is smaller than precision
1916 for (minstep = 0, i = 0; i < n; i++)
1926 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1928 if (stepsize < minstep)
1934 // Before taking any steps along the line, store the old position
1936 real *lastx = static_cast<real *>(last->s.x.data()[0]);
1937 real *lastf = static_cast<real *>(last->f.data()[0]);
1942 /* Take a step downhill.
1943 * In theory, we should find the actual minimum of the function in this
1944 * direction, somewhere along the line.
1945 * That is quite possible, but it turns out to take 5-10 function evaluations
1946 * for each line. However, we dont really need to find the exact minimum -
1947 * it is much better to start a new BFGS step in a modified direction as soon
1948 * as we are close to it. This will save a lot of energy evaluations.
1950 * In practice, we just try to take a single step.
1951 * If it worked (i.e. lowered the energy), we increase the stepsize but
1952 * continue straight to the next BFGS step without trying to find any minimum,
1953 * i.e. we change the search direction too. If the line was smooth, it is
1954 * likely we are in a smooth region, and then it makes sense to take longer
1955 * steps in the modified search direction too.
1957 * If it didn't work (higher energy), there must be a minimum somewhere between
1958 * the old position and the new one. Then we need to start by finding a lower
1959 * value before we change search direction. Since the energy was apparently
1960 * quite rough, we need to decrease the step size.
1962 * Due to the finite numerical accuracy, it turns out that it is a good idea
1963 * to accept a SMALL increase in energy, if the derivative is still downhill.
1964 * This leads to lower final energies in the tests I've done. / Erik
1967 // State "A" is the first position along the line.
1968 // reference position along line is initially zero
1971 // Check stepsize first. We do not allow displacements
1972 // larger than emstep.
1976 // Pick a new position C by adding stepsize to A.
1979 // Calculate what the largest change in any individual coordinate
1980 // would be (translation along line * gradient along line)
1982 for (i = 0; i < n; i++)
1985 if (delta > maxdelta)
1990 // If any displacement is larger than the stepsize limit, reduce the step
1991 if (maxdelta > inputrec->em_stepsize)
1996 while (maxdelta > inputrec->em_stepsize);
1998 // Take a trial step and move the coordinate array xc[] to position C
1999 real *xc = static_cast<real *>(sc->s.x.rvec_array()[0]);
2000 for (i = 0; i < n; i++)
2002 xc[i] = lastx[i] + c*s[i];
2006 // Calculate energy for the trial step in position C
2007 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2009 // Calc line gradient in position C
2010 real *fc = static_cast<real *>(sc->f.rvec_array()[0]);
2011 for (gpc = 0, i = 0; i < n; i++)
2013 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2015 /* Sum the gradient along the line across CPUs */
2018 gmx_sumd(1, &gpc, cr);
2021 // This is the max amount of increase in energy we tolerate.
2022 // By allowing VERY small changes (close to numerical precision) we
2023 // frequently find even better (lower) final energies.
2024 tmp = std::sqrt(GMX_REAL_EPS)*fabs(sa->epot);
2026 // Accept the step if the energy is lower in the new position C (compared to A),
2027 // or if it is not significantly higher and the line derivative is still negative.
2028 foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2029 // If true, great, we found a better energy. We no longer try to alter the
2030 // stepsize, but simply accept this new better position. The we select a new
2031 // search direction instead, which will be much more efficient than continuing
2032 // to take smaller steps along a line. Set fnorm based on the new C position,
2033 // which will be used to update the stepsize to 1/fnorm further down.
2035 // If false, the energy is NOT lower in point C, i.e. it will be the same
2036 // or higher than in point A. In this case it is pointless to move to point C,
2037 // so we will have to do more iterations along the same line to find a smaller
2038 // value in the interval [A=0.0,C].
2039 // Here, A is still 0.0, but that will change when we do a search in the interval
2040 // [0.0,C] below. That search we will do by interpolation or bisection rather
2041 // than with the stepsize, so no need to modify it. For the next search direction
2042 // it will be reset to 1/fnorm anyway.
2046 // OK, if we didn't find a lower value we will have to locate one now - there must
2047 // be one in the interval [a,c].
2048 // The same thing is valid here, though: Don't spend dozens of iterations to find
2049 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2050 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2051 // I also have a safeguard for potentially really pathological functions so we never
2052 // take more than 20 steps before we give up.
2053 // If we already found a lower value we just skip this step and continue to the update.
2058 // Select a new trial point B in the interval [A,C].
2059 // If the derivatives at points a & c have different sign we interpolate to zero,
2060 // otherwise just do a bisection since there might be multiple minima/maxima
2061 // inside the interval.
2062 if (gpa < 0 && gpc > 0)
2064 b = a + gpa*(a-c)/(gpc-gpa);
2071 /* safeguard if interpolation close to machine accuracy causes errors:
2072 * never go outside the interval
2074 if (b <= a || b >= c)
2079 // Take a trial step to point B
2080 real *xb = static_cast<real *>(sb->s.x.rvec_array()[0]);
2081 for (i = 0; i < n; i++)
2083 xb[i] = lastx[i] + b*s[i];
2087 // Calculate energy for the trial step in point B
2088 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2091 // Calculate gradient in point B
2092 real *fb = static_cast<real *>(sb->f.rvec_array()[0]);
2093 for (gpb = 0, i = 0; i < n; i++)
2095 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2098 /* Sum the gradient along the line across CPUs */
2101 gmx_sumd(1, &gpb, cr);
2104 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2105 // at the new point B, and rename the endpoints of this new interval A and C.
2108 /* Replace c endpoint with b */
2110 /* copy state b to c */
2115 /* Replace a endpoint with b */
2117 /* copy state b to a */
2122 * Stop search as soon as we find a value smaller than the endpoints,
2123 * or if the tolerance is below machine precision.
2124 * Never run more than 20 steps, no matter what.
2128 while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2130 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2132 /* OK. We couldn't find a significantly lower energy.
2133 * If ncorr==0 this was steepest descent, and then we give up.
2134 * If not, reset memory to restart as steepest descent before quitting.
2146 /* Search in gradient direction */
2147 for (i = 0; i < n; i++)
2149 dx[point][i] = ff[i];
2151 /* Reset stepsize */
2152 stepsize = 1.0/fnorm;
2157 /* Select min energy state of A & C, put the best in xx/ff/Epot
2159 if (sc->epot < sa->epot)
2181 /* Update the memory information, and calculate a new
2182 * approximation of the inverse hessian
2185 /* Have new data in Epot, xx, ff */
2186 if (ncorr < nmaxcorr)
2191 for (i = 0; i < n; i++)
2193 dg[point][i] = lastf[i]-ff[i];
2194 dx[point][i] *= step_taken;
2199 for (i = 0; i < n; i++)
2201 dgdg += dg[point][i]*dg[point][i];
2202 dgdx += dg[point][i]*dx[point][i];
2207 rho[point] = 1.0/dgdx;
2210 if (point >= nmaxcorr)
2216 for (i = 0; i < n; i++)
2223 /* Recursive update. First go back over the memory points */
2224 for (k = 0; k < ncorr; k++)
2233 for (i = 0; i < n; i++)
2235 sq += dx[cp][i]*p[i];
2238 alpha[cp] = rho[cp]*sq;
2240 for (i = 0; i < n; i++)
2242 p[i] -= alpha[cp]*dg[cp][i];
2246 for (i = 0; i < n; i++)
2251 /* And then go forward again */
2252 for (k = 0; k < ncorr; k++)
2255 for (i = 0; i < n; i++)
2257 yr += p[i]*dg[cp][i];
2261 beta = alpha[cp]-beta;
2263 for (i = 0; i < n; i++)
2265 p[i] += beta*dx[cp][i];
2275 for (i = 0; i < n; i++)
2279 dx[point][i] = p[i];
2287 /* Print it if necessary */
2290 if (mdrunOptions.verbose)
2292 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2293 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2294 step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2297 /* Store the new (lower) energies */
2298 matrix nullBox = {};
2299 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(step),
2300 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2301 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2303 do_log = do_per_step(step, inputrec->nstlog);
2304 do_ene = do_per_step(step, inputrec->nstenergy);
2306 imdSession->fillEnergyRecord(step, TRUE);
2310 energyOutput.printHeader(fplog, step, step);
2312 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2313 do_log ? fplog : nullptr, step, step,
2317 /* Send x and E to IMD client, if bIMD is TRUE. */
2318 if (imdSession->run(step, TRUE, state_global->box, state_global->x.rvec_array(), 0) && MASTER(cr))
2320 imdSession->sendPositionsAndEnergies();
2323 // Reset stepsize in we are doing more iterations
2326 /* Stop when the maximum force lies below tolerance.
2327 * If we have reached machine precision, converged is already set to true.
2329 converged = converged || (ems.fmax < inputrec->em_tol);
2331 } /* End of the loop */
2335 step--; /* we never took that last step in this case */
2338 if (ems.fmax > inputrec->em_tol)
2342 warn_step(fplog, inputrec->em_tol, ems.fmax,
2343 step-1 == number_steps, FALSE);
2348 /* If we printed energy and/or logfile last step (which was the last step)
2349 * we don't have to do it again, but otherwise print the final values.
2351 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2353 energyOutput.printHeader(fplog, step, step);
2355 if (!do_ene || !do_log) /* Write final energy file entries */
2357 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2358 !do_log ? fplog : nullptr, step, step,
2362 /* Print some stuff... */
2365 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2369 * For accurate normal mode calculation it is imperative that we
2370 * store the last conformation into the full precision binary trajectory.
2372 * However, we should only do it if we did NOT already write this step
2373 * above (which we did if do_x or do_f was true).
2375 do_x = !do_per_step(step, inputrec->nstxout);
2376 do_f = !do_per_step(step, inputrec->nstfout);
2377 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2378 top_global, inputrec, step,
2379 &ems, state_global, observablesHistory);
2383 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2384 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2385 number_steps, &ems, sqrtNumAtoms);
2386 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2387 number_steps, &ems, sqrtNumAtoms);
2389 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2392 finish_em(cr, outf, walltime_accounting, wcycle);
2394 /* To print the actual number of steps we needed somewhere */
2395 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2399 LegacySimulator::do_steep()
2401 const char *SD = "Steepest Descents";
2403 gmx_global_stat_t gstat;
2407 gmx_bool bDone, bAbort, do_x, do_f;
2412 int steps_accepted = 0;
2413 auto mdatoms = mdAtoms->mdatoms();
2415 GMX_LOG(mdlog.info).asParagraph().
2416 appendText("Note that activating steepest-descent energy minimization via the "
2417 "integrator .mdp option and the command gmx mdrun may "
2418 "be available in a different form in a future version of GROMACS, "
2419 "e.g. gmx minimize and an .mdp option.");
2421 /* Create 2 states on the stack and extract pointers that we will swap */
2422 em_state_t s0 {}, s1 {};
2423 em_state_t *s_min = &s0;
2424 em_state_t *s_try = &s1;
2426 /* Init em and store the local state in s_try */
2427 init_em(fplog, mdlog, SD, cr, inputrec, imdSession,
2429 state_global, top_global, s_try, &top,
2430 nrnb, fr, &graph, mdAtoms, &gstat,
2431 vsite, constr, nullptr);
2432 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2433 StartingBehavior::NewSimulation);
2434 gmx::EnergyOutput energyOutput(mdoutf_get_fp_ene(outf), top_global, inputrec, pull_work, nullptr, false, mdModulesNotifier);
2436 /* Print to log file */
2437 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2439 /* Set variables for stepsize (in nm). This is the largest
2440 * step that we are going to make in any direction.
2442 ustep = inputrec->em_stepsize;
2445 /* Max number of steps */
2446 nsteps = inputrec->nsteps;
2450 /* Print to the screen */
2451 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2455 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2457 EnergyEvaluator energyEvaluator {
2458 fplog, mdlog, cr, ms,
2460 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2461 vsite, constr, fcd, graph,
2462 mdAtoms, fr, runScheduleWork, enerd
2465 /**** HERE STARTS THE LOOP ****
2466 * count is the counter for the number of steps
2467 * bDone will be TRUE when the minimization has converged
2468 * bAbort will be TRUE when nsteps steps have been performed or when
2469 * the stepsize becomes smaller than is reasonable for machine precision
2474 while (!bDone && !bAbort)
2476 bAbort = (nsteps >= 0) && (count == nsteps);
2478 /* set new coordinates, except for first step */
2479 bool validStep = true;
2483 do_em_step(cr, inputrec, mdatoms,
2484 s_min, stepsize, &s_min->f, s_try,
2490 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2494 // Signal constraint error during stepping with energy=inf
2495 s_try->epot = std::numeric_limits<real>::infinity();
2500 energyOutput.printHeader(fplog, count, count);
2505 s_min->epot = s_try->epot;
2508 /* Print it if necessary */
2511 if (mdrunOptions.verbose)
2513 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2514 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2515 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2519 if ( (count == 0) || (s_try->epot < s_min->epot) )
2521 /* Store the new (lower) energies */
2522 matrix nullBox = {};
2523 energyOutput.addDataAtEnergyStep(false, false, static_cast<double>(count),
2524 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2525 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2527 imdSession->fillEnergyRecord(count, TRUE);
2529 const bool do_dr = do_per_step(steps_accepted, inputrec->nstdisreout);
2530 const bool do_or = do_per_step(steps_accepted, inputrec->nstorireout);
2531 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), TRUE,
2533 fplog, count, count,
2539 /* Now if the new energy is smaller than the previous...
2540 * or if this is the first step!
2541 * or if we did random steps!
2544 if ( (count == 0) || (s_try->epot < s_min->epot) )
2548 /* Test whether the convergence criterion is met... */
2549 bDone = (s_try->fmax < inputrec->em_tol);
2551 /* Copy the arrays for force, positions and energy */
2552 /* The 'Min' array always holds the coords and forces of the minimal
2554 swap_em_state(&s_min, &s_try);
2560 /* Write to trn, if necessary */
2561 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2562 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2563 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2564 top_global, inputrec, count,
2565 s_min, state_global, observablesHistory);
2569 /* If energy is not smaller make the step smaller... */
2572 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2574 /* Reload the old state */
2575 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec, imdSession,
2577 s_min, &top, mdAtoms, fr, vsite, constr,
2582 // If the force is very small after finishing minimization,
2583 // we risk dividing by zero when calculating the step size.
2584 // So we check first if the minimization has stopped before
2585 // trying to obtain a new step size.
2588 /* Determine new step */
2589 stepsize = ustep/s_min->fmax;
2592 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2594 if (count == nsteps || ustep < 1e-12)
2596 if (count == nsteps || ustep < 1e-6)
2601 warn_step(fplog, inputrec->em_tol, s_min->fmax,
2602 count == nsteps, constr != nullptr);
2607 /* Send IMD energies and positions, if bIMD is TRUE. */
2608 if (imdSession->run(count, TRUE, state_global->box,
2609 MASTER(cr) ? state_global->x.rvec_array() : nullptr,
2613 imdSession->sendPositionsAndEnergies();
2617 } /* End of the loop */
2619 /* Print some data... */
2622 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2624 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2625 top_global, inputrec, count,
2626 s_min, state_global, observablesHistory);
2630 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2632 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2633 s_min, sqrtNumAtoms);
2634 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2635 s_min, sqrtNumAtoms);
2638 finish_em(cr, outf, walltime_accounting, wcycle);
2640 /* To print the actual number of steps we needed somewhere */
2641 inputrec->nsteps = count;
2643 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2647 LegacySimulator::do_nm()
2649 const char *NM = "Normal Mode Analysis";
2652 gmx_global_stat_t gstat;
2657 gmx_bool bSparse; /* use sparse matrix storage format */
2659 gmx_sparsematrix_t * sparse_matrix = nullptr;
2660 real * full_matrix = nullptr;
2662 /* added with respect to mdrun */
2664 real der_range = 10.0*std::sqrt(GMX_REAL_EPS);
2666 bool bIsMaster = MASTER(cr);
2667 auto mdatoms = mdAtoms->mdatoms();
2669 GMX_LOG(mdlog.info).asParagraph().
2670 appendText("Note that activating normal-mode analysis via the integrator "
2671 ".mdp option and the command gmx mdrun may "
2672 "be available in a different form in a future version of GROMACS, "
2673 "e.g. gmx normal-modes.");
2675 if (constr != nullptr)
2677 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2680 gmx_shellfc_t *shellfc;
2682 em_state_t state_work {};
2684 /* Init em and store the local state in state_minimum */
2685 init_em(fplog, mdlog, NM, cr, inputrec, imdSession,
2687 state_global, top_global, &state_work, &top,
2688 nrnb, fr, &graph, mdAtoms, &gstat,
2689 vsite, constr, &shellfc);
2690 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, mdModulesNotifier, inputrec, top_global, nullptr, wcycle,
2691 StartingBehavior::NewSimulation);
2693 std::vector<int> atom_index = get_atom_index(top_global);
2694 std::vector<gmx::RVec> fneg(atom_index.size(), {0, 0, 0});
2695 snew(dfdx, atom_index.size());
2701 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2702 " which MIGHT not be accurate enough for normal mode analysis.\n"
2703 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2704 " are fairly modest even if you recompile in double precision.\n\n");
2708 /* Check if we can/should use sparse storage format.
2710 * Sparse format is only useful when the Hessian itself is sparse, which it
2711 * will be when we use a cutoff.
2712 * For small systems (n<1000) it is easier to always use full matrix format, though.
2714 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2716 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2719 else if (atom_index.size() < 1000)
2721 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2727 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2731 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2732 sz = DIM*atom_index.size();
2734 fprintf(stderr, "Allocating Hessian memory...\n\n");
2738 sparse_matrix = gmx_sparsematrix_init(sz);
2739 sparse_matrix->compressed_symmetric = TRUE;
2743 snew(full_matrix, sz*sz);
2746 /* Write start time and temperature */
2747 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2749 /* fudge nr of steps to nr of atoms */
2750 inputrec->nsteps = atom_index.size()*2;
2754 fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2755 *(top_global->name), inputrec->nsteps);
2758 nnodes = cr->nnodes;
2760 /* Make evaluate_energy do a single node force calculation */
2762 EnergyEvaluator energyEvaluator {
2763 fplog, mdlog, cr, ms,
2765 inputrec, imdSession, pull_work, nrnb, wcycle, gstat,
2766 vsite, constr, fcd, graph,
2767 mdAtoms, fr, runScheduleWork, enerd
2769 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2770 cr->nnodes = nnodes;
2772 /* if forces are not small, warn user */
2773 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2775 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2776 if (state_work.fmax > 1.0e-3)
2778 GMX_LOG(mdlog.warning).appendText(
2779 "The force is probably not small enough to "
2780 "ensure that you are at a minimum.\n"
2781 "Be aware that negative eigenvalues may occur\n"
2782 "when the resulting matrix is diagonalized.");
2785 /***********************************************************
2787 * Loop over all pairs in matrix
2789 * do_force called twice. Once with positive and
2790 * once with negative displacement
2792 ************************************************************/
2794 /* Steps are divided one by one over the nodes */
2796 auto state_work_x = makeArrayRef(state_work.s.x);
2797 auto state_work_f = makeArrayRef(state_work.f);
2798 for (index aid = cr->nodeid; aid < ssize(atom_index); aid += nnodes)
2800 size_t atom = atom_index[aid];
2801 for (size_t d = 0; d < DIM; d++)
2804 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2807 x_min = state_work_x[atom][d];
2809 for (unsigned int dx = 0; (dx < 2); dx++)
2813 state_work_x[atom][d] = x_min - der_range;
2817 state_work_x[atom][d] = x_min + der_range;
2820 /* Make evaluate_energy do a single node force calculation */
2824 /* Now is the time to relax the shells */
2825 relax_shell_flexcon(fplog,
2828 mdrunOptions.verbose,
2840 state_work.s.natoms,
2841 state_work.s.x.arrayRefWithPadding(),
2842 state_work.s.v.arrayRefWithPadding(),
2844 state_work.s.lambda,
2846 state_work.f.arrayRefWithPadding(),
2858 DDBalanceRegionHandler(nullptr));
2864 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid*2+dx, FALSE);
2867 cr->nnodes = nnodes;
2871 std::copy(state_work_f.begin(), state_work_f.begin()+atom_index.size(), fneg.begin());
2875 /* x is restored to original */
2876 state_work_x[atom][d] = x_min;
2878 for (size_t j = 0; j < atom_index.size(); j++)
2880 for (size_t k = 0; (k < DIM); k++)
2883 -(state_work_f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2890 #define mpi_type GMX_MPI_REAL
2891 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2892 cr->nodeid, cr->mpi_comm_mygroup);
2897 for (index node = 0; (node < nnodes && aid+node < ssize(atom_index)); node++)
2903 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2904 cr->mpi_comm_mygroup, &stat);
2909 row = (aid + node)*DIM + d;
2911 for (size_t j = 0; j < atom_index.size(); j++)
2913 for (size_t k = 0; k < DIM; k++)
2919 if (col >= row && dfdx[j][k] != 0.0)
2921 gmx_sparsematrix_increment_value(sparse_matrix,
2922 row, col, dfdx[j][k]);
2927 full_matrix[row*sz+col] = dfdx[j][k];
2934 if (mdrunOptions.verbose && fplog)
2939 /* write progress */
2940 if (bIsMaster && mdrunOptions.verbose)
2942 fprintf(stderr, "\rFinished step %d out of %td",
2943 std::min<int>(atom+nnodes, atom_index.size()),
2951 fprintf(stderr, "\n\nWriting Hessian...\n");
2952 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2955 finish_em(cr, outf, walltime_accounting, wcycle);
2957 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);