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/domdec.h"
59 #include "gromacs/domdec/domdec_struct.h"
60 #include "gromacs/domdec/partition.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fileio/confio.h"
63 #include "gromacs/fileio/mtxio.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/imd/imd.h"
67 #include "gromacs/linearalgebra/sparsematrix.h"
68 #include "gromacs/listed_forces/manage_threading.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/vec.h"
71 #include "gromacs/mdlib/constr.h"
72 #include "gromacs/mdlib/force.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/gmx_omp_nthreads.h"
75 #include "gromacs/mdlib/md_support.h"
76 #include "gromacs/mdlib/mdatoms.h"
77 #include "gromacs/mdlib/mdebin.h"
78 #include "gromacs/mdlib/mdrun.h"
79 #include "gromacs/mdlib/mdsetup.h"
80 #include "gromacs/mdlib/ns.h"
81 #include "gromacs/mdlib/shellfc.h"
82 #include "gromacs/mdlib/sim_util.h"
83 #include "gromacs/mdlib/tgroup.h"
84 #include "gromacs/mdlib/trajectory_writing.h"
85 #include "gromacs/mdlib/update.h"
86 #include "gromacs/mdlib/vsite.h"
87 #include "gromacs/mdtypes/commrec.h"
88 #include "gromacs/mdtypes/inputrec.h"
89 #include "gromacs/mdtypes/md_enums.h"
90 #include "gromacs/mdtypes/state.h"
91 #include "gromacs/pbcutil/mshift.h"
92 #include "gromacs/pbcutil/pbc.h"
93 #include "gromacs/timing/wallcycle.h"
94 #include "gromacs/timing/walltime_accounting.h"
95 #include "gromacs/topology/mtop_util.h"
96 #include "gromacs/topology/topology.h"
97 #include "gromacs/utility/cstringutil.h"
98 #include "gromacs/utility/exceptions.h"
99 #include "gromacs/utility/fatalerror.h"
100 #include "gromacs/utility/logger.h"
101 #include "gromacs/utility/smalloc.h"
103 #include "integrator.h"
105 //! Utility structure for manipulating states during EM
107 //! Copy of the global state
110 PaddedVector<gmx::RVec> f;
113 //! Norm of the force
121 //! Print the EM starting conditions
122 static void print_em_start(FILE *fplog,
124 gmx_walltime_accounting_t walltime_accounting,
125 gmx_wallcycle_t wcycle,
128 walltime_accounting_start_time(walltime_accounting);
129 wallcycle_start(wcycle, ewcRUN);
130 print_start(fplog, cr, walltime_accounting, name);
133 //! Stop counting time for EM
134 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
135 gmx_wallcycle_t wcycle)
137 wallcycle_stop(wcycle, ewcRUN);
139 walltime_accounting_end_time(walltime_accounting);
142 //! Printing a log file and console header
143 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
146 fprintf(out, "%s:\n", minimizer);
147 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
148 fprintf(out, " Number of steps = %12d\n", nsteps);
151 //! Print warning message
152 static void warn_step(FILE *fp,
158 constexpr bool realIsDouble = GMX_DOUBLE;
161 if (!std::isfinite(fmax))
164 "\nEnergy minimization has stopped because the force "
165 "on at least one atom is not finite. This usually means "
166 "atoms are overlapping. Modify the input coordinates to "
167 "remove atom overlap or use soft-core potentials with "
168 "the free energy code to avoid infinite forces.\n%s",
170 "You could also be lucky that switching to double precision "
171 "is sufficient to obtain finite forces.\n" :
177 "\nEnergy minimization reached the maximum number "
178 "of steps before the forces reached the requested "
179 "precision Fmax < %g.\n", ftol);
184 "\nEnergy minimization has stopped, but the forces have "
185 "not converged to the requested precision Fmax < %g (which "
186 "may not be possible for your system). It stopped "
187 "because the algorithm tried to make a new step whose size "
188 "was too small, or there was no change in the energy since "
189 "last step. Either way, we regard the minimization as "
190 "converged to within the available machine precision, "
191 "given your starting configuration and EM parameters.\n%s%s",
194 "\nDouble precision normally gives you higher accuracy, but "
195 "this is often not needed for preparing to run molecular "
199 "You might need to increase your constraint accuracy, or turn\n"
200 "off constraints altogether (set constraints = none in mdp file)\n" :
204 fputs(wrap_lines(buffer, 78, 0, FALSE), stderr);
205 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
208 //! Print message about convergence of the EM
209 static void print_converged(FILE *fp, const char *alg, real ftol,
210 int64_t count, gmx_bool bDone, int64_t nsteps,
211 const em_state_t *ems, double sqrtNumAtoms)
213 char buf[STEPSTRSIZE];
217 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
218 alg, ftol, gmx_step_str(count, buf));
220 else if (count < nsteps)
222 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
223 "but did not reach the requested Fmax < %g.\n",
224 alg, gmx_step_str(count, buf), ftol);
228 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
229 alg, ftol, gmx_step_str(count, buf));
233 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
234 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
235 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm/sqrtNumAtoms);
237 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
238 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
239 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm/sqrtNumAtoms);
243 //! Compute the norm and max of the force array in parallel
244 static void get_f_norm_max(const t_commrec *cr,
245 t_grpopts *opts, t_mdatoms *mdatoms, const rvec *f,
246 real *fnorm, real *fmax, int *a_fmax)
250 int la_max, a_max, start, end, i, m, gf;
252 /* This routine finds the largest force and returns it.
253 * On parallel machines the global max is taken.
259 end = mdatoms->homenr;
260 if (mdatoms->cFREEZE)
262 for (i = start; i < end; i++)
264 gf = mdatoms->cFREEZE[i];
266 for (m = 0; m < DIM; m++)
268 if (!opts->nFreeze[gf][m])
270 fam += gmx::square(f[i][m]);
283 for (i = start; i < end; i++)
295 if (la_max >= 0 && DOMAINDECOMP(cr))
297 a_max = cr->dd->globalAtomIndices[la_max];
305 snew(sum, 2*cr->nnodes+1);
306 sum[2*cr->nodeid] = fmax2;
307 sum[2*cr->nodeid+1] = a_max;
308 sum[2*cr->nnodes] = fnorm2;
309 gmx_sumd(2*cr->nnodes+1, sum, cr);
310 fnorm2 = sum[2*cr->nnodes];
311 /* Determine the global maximum */
312 for (i = 0; i < cr->nnodes; i++)
314 if (sum[2*i] > fmax2)
317 a_max = gmx::roundToInt(sum[2*i+1]);
325 *fnorm = sqrt(fnorm2);
337 //! Compute the norm of the force
338 static void get_state_f_norm_max(const t_commrec *cr,
339 t_grpopts *opts, t_mdatoms *mdatoms,
342 get_f_norm_max(cr, opts, mdatoms, ems->f.rvec_array(),
343 &ems->fnorm, &ems->fmax, &ems->a_fmax);
346 //! Initialize the energy minimization
347 static void init_em(FILE *fplog,
348 const gmx::MDLogger &mdlog,
351 const gmx_multisim_t *ms,
353 const MdrunOptions &mdrunOptions,
354 t_state *state_global, gmx_mtop_t *top_global,
355 em_state_t *ems, gmx_localtop_t *top,
358 t_graph **graph, gmx::MDAtoms *mdAtoms, gmx_global_stat_t *gstat,
359 gmx_vsite_t *vsite, gmx::Constraints *constr, gmx_shellfc_t **shellfc,
360 int nfile, const t_filenm fnm[])
366 fprintf(fplog, "Initiating %s\n", title);
371 state_global->ngtc = 0;
373 initialize_lambdas(fplog, *ir, MASTER(cr), &(state_global->fep_state), state_global->lambda, nullptr);
377 /* Interactive molecular dynamics */
378 init_IMD(ir, cr, ms, top_global, fplog, 1,
379 MASTER(cr) ? state_global->x.rvec_array() : nullptr,
380 nfile, fnm, nullptr, mdrunOptions);
384 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
386 *shellfc = init_shell_flexcon(stdout,
388 constr ? constr->numFlexibleConstraints() : 0,
394 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
396 /* With energy minimization, shells and flexible constraints are
397 * automatically minimized when treated like normal DOFS.
399 if (shellfc != nullptr)
405 auto mdatoms = mdAtoms->mdatoms();
406 if (DOMAINDECOMP(cr))
408 top->useInDomainDecomp_ = true;
409 dd_init_local_top(*top_global, top);
411 dd_init_local_state(cr->dd, state_global, &ems->s);
413 /* Distribute the charge groups over the nodes from the master node */
414 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
415 state_global, *top_global, ir,
416 &ems->s, &ems->f, mdAtoms, top,
418 nrnb, nullptr, FALSE);
419 dd_store_state(cr->dd, &ems->s);
425 state_change_natoms(state_global, state_global->natoms);
426 /* Just copy the state */
427 ems->s = *state_global;
428 state_change_natoms(&ems->s, ems->s.natoms);
429 ems->f.resizeWithPadding(ems->s.natoms);
431 mdAlgorithmsSetupAtomData(cr, ir, *top_global, top, fr,
433 constr, vsite, shellfc ? *shellfc : nullptr);
437 set_vsite_top(vsite, top, mdatoms);
441 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
445 // TODO how should this cross-module support dependency be managed?
446 if (ir->eConstrAlg == econtSHAKE &&
447 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
449 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
450 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
453 if (!ir->bContinuation)
455 /* Constrain the starting coordinates */
457 constr->apply(TRUE, TRUE,
459 ems->s.x.rvec_array(),
460 ems->s.x.rvec_array(),
463 ems->s.lambda[efptFEP], &dvdl_constr,
464 nullptr, nullptr, gmx::ConstraintVariable::Positions);
470 *gstat = global_stat_init(ir);
477 calc_shifts(ems->s.box, fr->shift_vec);
480 //! Finalize the minimization
481 static void finish_em(const t_commrec *cr, gmx_mdoutf_t outf,
482 gmx_walltime_accounting_t walltime_accounting,
483 gmx_wallcycle_t wcycle)
485 if (!thisRankHasDuty(cr, DUTY_PME))
487 /* Tell the PME only node to finish */
488 gmx_pme_send_finish(cr);
493 em_time_end(walltime_accounting, wcycle);
496 //! Swap two different EM states during minimization
497 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
506 //! Save the EM trajectory
507 static void write_em_traj(FILE *fplog, const t_commrec *cr,
509 gmx_bool bX, gmx_bool bF, const char *confout,
510 gmx_mtop_t *top_global,
511 t_inputrec *ir, int64_t step,
513 t_state *state_global,
514 ObservablesHistory *observablesHistory)
520 mdof_flags |= MDOF_X;
524 mdof_flags |= MDOF_F;
527 /* If we want IMD output, set appropriate MDOF flag */
530 mdof_flags |= MDOF_IMD;
533 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
534 top_global, step, static_cast<double>(step),
535 &state->s, state_global, observablesHistory,
538 if (confout != nullptr)
540 if (DOMAINDECOMP(cr))
542 /* If bX=true, x was collected to state_global in the call above */
545 gmx::ArrayRef<gmx::RVec> globalXRef = MASTER(cr) ? makeArrayRef(state_global->x) : gmx::EmptyArrayRef();
546 dd_collect_vec(cr->dd, &state->s, makeArrayRef(state->s.x), globalXRef);
551 /* Copy the local state pointer */
552 state_global = &state->s;
557 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
559 /* Make molecules whole only for confout writing */
560 do_pbc_mtop(fplog, ir->ePBC, state->s.box, top_global,
561 state_global->x.rvec_array());
564 write_sto_conf_mtop(confout,
565 *top_global->name, top_global,
566 state_global->x.rvec_array(), nullptr, ir->ePBC, state->s.box);
571 //! \brief Do one minimization step
573 // \returns true when the step succeeded, false when a constraint error occurred
574 static bool do_em_step(const t_commrec *cr,
575 t_inputrec *ir, t_mdatoms *md,
576 em_state_t *ems1, real a, const PaddedVector<gmx::RVec> *force,
578 gmx::Constraints *constr,
585 int nthreads gmx_unused;
587 bool validStep = true;
592 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
594 gmx_incons("state mismatch in do_em_step");
597 s2->flags = s1->flags;
599 if (s2->natoms != s1->natoms)
601 state_change_natoms(s2, s1->natoms);
602 ems2->f.resizeWithPadding(s2->natoms);
604 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
606 s2->cg_gl.resize(s1->cg_gl.size());
609 copy_mat(s1->box, s2->box);
610 /* Copy free energy state */
611 s2->lambda = s1->lambda;
612 copy_mat(s1->box, s2->box);
617 nthreads = gmx_omp_nthreads_get(emntUpdate);
618 #pragma omp parallel num_threads(nthreads)
620 const rvec *x1 = s1->x.rvec_array();
621 rvec *x2 = s2->x.rvec_array();
622 const rvec *f = force->rvec_array();
625 #pragma omp for schedule(static) nowait
626 for (int i = start; i < end; i++)
634 for (int m = 0; m < DIM; m++)
636 if (ir->opts.nFreeze[gf][m])
642 x2[i][m] = x1[i][m] + a*f[i][m];
646 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
649 if (s2->flags & (1<<estCGP))
651 /* Copy the CG p vector */
652 const rvec *p1 = s1->cg_p.rvec_array();
653 rvec *p2 = s2->cg_p.rvec_array();
654 #pragma omp for schedule(static) nowait
655 for (int i = start; i < end; i++)
657 // Trivial OpenMP block that does not throw
658 copy_rvec(p1[i], p2[i]);
662 if (DOMAINDECOMP(cr))
664 s2->ddp_count = s1->ddp_count;
666 /* OpenMP does not supported unsigned loop variables */
667 #pragma omp for schedule(static) nowait
668 for (int i = 0; i < gmx::ssize(s2->cg_gl); i++)
670 s2->cg_gl[i] = s1->cg_gl[i];
672 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
680 constr->apply(TRUE, TRUE,
682 s1->x.rvec_array(), s2->x.rvec_array(),
684 s2->lambda[efptBONDED], &dvdl_constr,
685 nullptr, nullptr, gmx::ConstraintVariable::Positions);
689 /* This global reduction will affect performance at high
690 * parallelization, but we can not really avoid it.
691 * But usually EM is not run at high parallelization.
693 int reductionBuffer = static_cast<int>(!validStep);
694 gmx_sumi(1, &reductionBuffer, cr);
695 validStep = (reductionBuffer == 0);
698 // We should move this check to the different minimizers
699 if (!validStep && ir->eI != eiSteep)
701 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
702 EI(ir->eI), EI(eiSteep), EI(ir->eI));
709 //! Prepare EM for using domain decomposition parallellization
710 static void em_dd_partition_system(FILE *fplog,
711 const gmx::MDLogger &mdlog,
712 int step, const t_commrec *cr,
713 gmx_mtop_t *top_global, t_inputrec *ir,
714 em_state_t *ems, gmx_localtop_t *top,
715 gmx::MDAtoms *mdAtoms, t_forcerec *fr,
716 gmx_vsite_t *vsite, gmx::Constraints *constr,
717 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
719 /* Repartition the domain decomposition */
720 dd_partition_system(fplog, mdlog, step, cr, FALSE, 1,
721 nullptr, *top_global, ir,
723 mdAtoms, top, fr, vsite, constr,
724 nrnb, wcycle, FALSE);
725 dd_store_state(cr->dd, &ems->s);
731 /*! \brief Class to handle the work of setting and doing an energy evaluation.
733 * This class is a mere aggregate of parameters to pass to evaluate an
734 * energy, so that future changes to names and types of them consume
735 * less time when refactoring other code.
737 * Aggregate initialization is used, for which the chief risk is that
738 * if a member is added at the end and not all initializer lists are
739 * updated, then the member will be value initialized, which will
740 * typically mean initialization to zero.
742 * We only want to construct one of these with an initializer list, so
743 * we explicitly delete the default constructor. */
744 class EnergyEvaluator
747 //! We only intend to construct such objects with an initializer list.
748 EnergyEvaluator() = delete;
749 /*! \brief Evaluates an energy on the state in \c ems.
751 * \todo In practice, the same objects mu_tot, vir, and pres
752 * are always passed to this function, so we would rather have
753 * them as data members. However, their C-array types are
754 * unsuited for aggregate initialization. When the types
755 * improve, the call signature of this method can be reduced.
757 void run(em_state_t *ems, rvec mu_tot,
758 tensor vir, tensor pres,
759 int64_t count, gmx_bool bFirst);
760 //! Handles logging (deprecated).
763 const gmx::MDLogger &mdlog;
764 //! Handles communication.
766 //! Coordinates multi-simulations.
767 const gmx_multisim_t *ms;
768 //! Holds the simulation topology.
769 gmx_mtop_t *top_global;
770 //! Holds the domain topology.
772 //! User input options.
773 t_inputrec *inputrec;
774 //! Manages flop accounting.
776 //! Manages wall cycle accounting.
777 gmx_wallcycle_t wcycle;
778 //! Coordinates global reduction.
779 gmx_global_stat_t gstat;
780 //! Handles virtual sites.
782 //! Handles constraints.
783 gmx::Constraints *constr;
784 //! Handles strange things.
786 //! Molecular graph for SHAKE.
788 //! Per-atom data for this domain.
789 gmx::MDAtoms *mdAtoms;
790 //! Handles how to calculate the forces.
792 //! Schedule of force-calculation work each step for this task.
793 gmx::PpForceWorkload *ppForceWorkload;
794 //! Stores the computed energies.
795 gmx_enerdata_t *enerd;
799 EnergyEvaluator::run(em_state_t *ems, rvec mu_tot,
800 tensor vir, tensor pres,
801 int64_t count, gmx_bool bFirst)
805 tensor force_vir, shake_vir, ekin;
806 real dvdl_constr, prescorr, enercorr, dvdlcorr;
809 /* Set the time to the initial time, the time does not change during EM */
810 t = inputrec->init_t;
813 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
815 /* This is the first state or an old state used before the last ns */
821 if (inputrec->nstlist > 0)
829 construct_vsites(vsite, ems->s.x.rvec_array(), 1, nullptr,
830 top->idef.iparams, top->idef.il,
831 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
834 if (DOMAINDECOMP(cr) && bNS)
836 /* Repartition the domain decomposition */
837 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec,
838 ems, top, mdAtoms, fr, vsite, constr,
842 /* Calc force & energy on new trial position */
843 /* do_force always puts the charge groups in the box and shifts again
844 * We do not unshift, so molecules are always whole in congrad.c
846 do_force(fplog, cr, ms, inputrec, nullptr, nullptr,
847 count, nrnb, wcycle, top, &top_global->groups,
848 ems->s.box, ems->s.x.arrayRefWithPadding(), &ems->s.hist,
849 ems->f.arrayRefWithPadding(), force_vir, mdAtoms->mdatoms(), enerd, fcd,
850 ems->s.lambda, graph, fr, ppForceWorkload, vsite, mu_tot, t, nullptr,
851 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
852 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
853 (bNS ? GMX_FORCE_NS : 0),
855 DdOpenBalanceRegionBeforeForceComputation::yes :
856 DdOpenBalanceRegionBeforeForceComputation::no,
858 DdCloseBalanceRegionAfterForceComputation::yes :
859 DdCloseBalanceRegionAfterForceComputation::no);
861 /* Clear the unused shake virial and pressure */
862 clear_mat(shake_vir);
865 /* Communicate stuff when parallel */
866 if (PAR(cr) && inputrec->eI != eiNM)
868 wallcycle_start(wcycle, ewcMoveE);
870 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
871 inputrec, nullptr, nullptr, nullptr, 1, &terminate,
877 wallcycle_stop(wcycle, ewcMoveE);
880 /* Calculate long range corrections to pressure and energy */
881 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
882 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
883 enerd->term[F_DISPCORR] = enercorr;
884 enerd->term[F_EPOT] += enercorr;
885 enerd->term[F_PRES] += prescorr;
886 enerd->term[F_DVDL] += dvdlcorr;
888 ems->epot = enerd->term[F_EPOT];
892 /* Project out the constraint components of the force */
894 rvec *f_rvec = ems->f.rvec_array();
895 constr->apply(FALSE, FALSE,
897 ems->s.x.rvec_array(), f_rvec, f_rvec,
899 ems->s.lambda[efptBONDED], &dvdl_constr,
900 nullptr, &shake_vir, gmx::ConstraintVariable::ForceDispl);
901 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
902 m_add(force_vir, shake_vir, vir);
906 copy_mat(force_vir, vir);
910 enerd->term[F_PRES] =
911 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
913 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
915 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
917 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
923 //! Parallel utility summing energies and forces
924 static double reorder_partsum(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
925 gmx_mtop_t *top_global,
926 em_state_t *s_min, em_state_t *s_b)
929 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
931 unsigned char *grpnrFREEZE;
935 fprintf(debug, "Doing reorder_partsum\n");
938 const rvec *fm = s_min->f.rvec_array();
939 const rvec *fb = s_b->f.rvec_array();
941 cgs_gl = dd_charge_groups_global(cr->dd);
942 index = cgs_gl->index;
944 /* Collect fm in a global vector fmg.
945 * This conflicts with the spirit of domain decomposition,
946 * but to fully optimize this a much more complicated algorithm is required.
949 snew(fmg, top_global->natoms);
951 ncg = s_min->s.cg_gl.size();
952 cg_gl = s_min->s.cg_gl.data();
954 for (c = 0; c < ncg; c++)
959 for (a = a0; a < a1; a++)
961 copy_rvec(fm[i], fmg[a]);
965 gmx_sum(top_global->natoms*3, fmg[0], cr);
967 /* Now we will determine the part of the sum for the cgs in state s_b */
968 ncg = s_b->s.cg_gl.size();
969 cg_gl = s_b->s.cg_gl.data();
973 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
974 for (c = 0; c < ncg; c++)
979 for (a = a0; a < a1; a++)
981 if (mdatoms->cFREEZE && grpnrFREEZE)
985 for (m = 0; m < DIM; m++)
987 if (!opts->nFreeze[gf][m])
989 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
1001 //! Print some stuff, like beta, whatever that means.
1002 static real pr_beta(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
1003 gmx_mtop_t *top_global,
1004 em_state_t *s_min, em_state_t *s_b)
1008 /* This is just the classical Polak-Ribiere calculation of beta;
1009 * it looks a bit complicated since we take freeze groups into account,
1010 * and might have to sum it in parallel runs.
1013 if (!DOMAINDECOMP(cr) ||
1014 (s_min->s.ddp_count == cr->dd->ddp_count &&
1015 s_b->s.ddp_count == cr->dd->ddp_count))
1017 const rvec *fm = s_min->f.rvec_array();
1018 const rvec *fb = s_b->f.rvec_array();
1021 /* This part of code can be incorrect with DD,
1022 * since the atom ordering in s_b and s_min might differ.
1024 for (int i = 0; i < mdatoms->homenr; i++)
1026 if (mdatoms->cFREEZE)
1028 gf = mdatoms->cFREEZE[i];
1030 for (int m = 0; m < DIM; m++)
1032 if (!opts->nFreeze[gf][m])
1034 sum += (fb[i][m] - fm[i][m])*fb[i][m];
1041 /* We need to reorder cgs while summing */
1042 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
1046 gmx_sumd(1, &sum, cr);
1049 return sum/gmx::square(s_min->fnorm);
1058 const char *CG = "Polak-Ribiere Conjugate Gradients";
1061 gmx_enerdata_t *enerd;
1062 gmx_global_stat_t gstat;
1064 double tmp, minstep;
1066 real a, b, c, beta = 0.0;
1069 gmx_bool converged, foundlower;
1071 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1073 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1074 int m, step, nminstep;
1075 auto mdatoms = mdAtoms->mdatoms();
1077 GMX_LOG(mdlog.info).asParagraph().
1078 appendText("Note that activating conjugate gradient energy minimization via the "
1079 "integrator .mdp option and the command gmx mdrun may "
1080 "be available in a different form in a future version of GROMACS, "
1081 "e.g. gmx minimize and an .mdp option.");
1087 // In CG, the state is extended with a search direction
1088 state_global->flags |= (1<<estCGP);
1090 // Ensure the extra per-atom state array gets allocated
1091 state_change_natoms(state_global, state_global->natoms);
1093 // Initialize the search direction to zero
1094 for (RVec &cg_p : state_global->cg_p)
1100 /* Create 4 states on the stack and extract pointers that we will swap */
1101 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1102 em_state_t *s_min = &s0;
1103 em_state_t *s_a = &s1;
1104 em_state_t *s_b = &s2;
1105 em_state_t *s_c = &s3;
1107 /* Init em and store the local state in s_min */
1108 init_em(fplog, mdlog, CG, cr, ms, inputrec, mdrunOptions,
1109 state_global, top_global, s_min, &top,
1110 nrnb, fr, &graph, mdAtoms, &gstat,
1111 vsite, constr, nullptr,
1113 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, inputrec, top_global, nullptr, wcycle);
1115 init_enerdata(top_global->groups.grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
1116 t_mdebin *mdebin = init_mdebin(mdoutf_get_fp_ene(outf), top_global, inputrec, nullptr);
1118 /* Print to log file */
1119 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1121 /* Max number of steps */
1122 number_steps = inputrec->nsteps;
1126 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1130 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1133 EnergyEvaluator energyEvaluator {
1134 fplog, mdlog, cr, ms,
1136 inputrec, nrnb, wcycle, gstat,
1137 vsite, constr, fcd, graph,
1138 mdAtoms, fr, ppForceWorkload, enerd
1140 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1141 /* do_force always puts the charge groups in the box and shifts again
1142 * We do not unshift, so molecules are always whole in congrad.c
1144 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1148 /* Copy stuff to the energy bin for easy printing etc. */
1149 matrix nullBox = {};
1150 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(step),
1151 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1152 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1154 print_ebin_header(fplog, step, step);
1155 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1156 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1159 /* Estimate/guess the initial stepsize */
1160 stepsize = inputrec->em_stepsize/s_min->fnorm;
1164 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1165 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1166 s_min->fmax, s_min->a_fmax+1);
1167 fprintf(stderr, " F-Norm = %12.5e\n",
1168 s_min->fnorm/sqrtNumAtoms);
1169 fprintf(stderr, "\n");
1170 /* and copy to the log file too... */
1171 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1172 s_min->fmax, s_min->a_fmax+1);
1173 fprintf(fplog, " F-Norm = %12.5e\n",
1174 s_min->fnorm/sqrtNumAtoms);
1175 fprintf(fplog, "\n");
1177 /* Start the loop over CG steps.
1178 * Each successful step is counted, and we continue until
1179 * we either converge or reach the max number of steps.
1182 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1185 /* start taking steps in a new direction
1186 * First time we enter the routine, beta=0, and the direction is
1187 * simply the negative gradient.
1190 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1191 rvec *pm = s_min->s.cg_p.rvec_array();
1192 const rvec *sfm = s_min->f.rvec_array();
1195 for (int i = 0; i < mdatoms->homenr; i++)
1197 if (mdatoms->cFREEZE)
1199 gf = mdatoms->cFREEZE[i];
1201 for (m = 0; m < DIM; m++)
1203 if (!inputrec->opts.nFreeze[gf][m])
1205 pm[i][m] = sfm[i][m] + beta*pm[i][m];
1206 gpa -= pm[i][m]*sfm[i][m];
1207 /* f is negative gradient, thus the sign */
1216 /* Sum the gradient along the line across CPUs */
1219 gmx_sumd(1, &gpa, cr);
1222 /* Calculate the norm of the search vector */
1223 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1225 /* Just in case stepsize reaches zero due to numerical precision... */
1228 stepsize = inputrec->em_stepsize/pnorm;
1232 * Double check the value of the derivative in the search direction.
1233 * If it is positive it must be due to the old information in the
1234 * CG formula, so just remove that and start over with beta=0.
1235 * This corresponds to a steepest descent step.
1240 step--; /* Don't count this step since we are restarting */
1241 continue; /* Go back to the beginning of the big for-loop */
1244 /* Calculate minimum allowed stepsize, before the average (norm)
1245 * relative change in coordinate is smaller than precision
1248 auto s_min_x = makeArrayRef(s_min->s.x);
1249 for (int i = 0; i < mdatoms->homenr; i++)
1251 for (m = 0; m < DIM; m++)
1253 tmp = fabs(s_min_x[i][m]);
1262 /* Add up from all CPUs */
1265 gmx_sumd(1, &minstep, cr);
1268 minstep = GMX_REAL_EPS/sqrt(minstep/(3*top_global->natoms));
1270 if (stepsize < minstep)
1276 /* Write coordinates if necessary */
1277 do_x = do_per_step(step, inputrec->nstxout);
1278 do_f = do_per_step(step, inputrec->nstfout);
1280 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1281 top_global, inputrec, step,
1282 s_min, state_global, observablesHistory);
1284 /* Take a step downhill.
1285 * In theory, we should minimize the function along this direction.
1286 * That is quite possible, but it turns out to take 5-10 function evaluations
1287 * for each line. However, we dont really need to find the exact minimum -
1288 * it is much better to start a new CG step in a modified direction as soon
1289 * as we are close to it. This will save a lot of energy evaluations.
1291 * In practice, we just try to take a single step.
1292 * If it worked (i.e. lowered the energy), we increase the stepsize but
1293 * the continue straight to the next CG step without trying to find any minimum.
1294 * If it didn't work (higher energy), there must be a minimum somewhere between
1295 * the old position and the new one.
1297 * Due to the finite numerical accuracy, it turns out that it is a good idea
1298 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1299 * This leads to lower final energies in the tests I've done. / Erik
1301 s_a->epot = s_min->epot;
1303 c = a + stepsize; /* reference position along line is zero */
1305 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1307 em_dd_partition_system(fplog, mdlog, step, cr, top_global, inputrec,
1308 s_min, &top, mdAtoms, fr, vsite, constr,
1312 /* Take a trial step (new coords in s_c) */
1313 do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c,
1317 /* Calculate energy for the trial step */
1318 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1320 /* Calc derivative along line */
1321 const rvec *pc = s_c->s.cg_p.rvec_array();
1322 const rvec *sfc = s_c->f.rvec_array();
1324 for (int i = 0; i < mdatoms->homenr; i++)
1326 for (m = 0; m < DIM; m++)
1328 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1331 /* Sum the gradient along the line across CPUs */
1334 gmx_sumd(1, &gpc, cr);
1337 /* This is the max amount of increase in energy we tolerate */
1338 tmp = std::sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1340 /* Accept the step if the energy is lower, or if it is not significantly higher
1341 * and the line derivative is still negative.
1343 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1346 /* Great, we found a better energy. Increase step for next iteration
1347 * if we are still going down, decrease it otherwise
1351 stepsize *= 1.618034; /* The golden section */
1355 stepsize *= 0.618034; /* 1/golden section */
1360 /* New energy is the same or higher. We will have to do some work
1361 * to find a smaller value in the interval. Take smaller step next time!
1364 stepsize *= 0.618034;
1370 /* OK, if we didn't find a lower value we will have to locate one now - there must
1371 * be one in the interval [a=0,c].
1372 * The same thing is valid here, though: Don't spend dozens of iterations to find
1373 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1374 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1376 * I also have a safeguard for potentially really pathological functions so we never
1377 * take more than 20 steps before we give up ...
1379 * If we already found a lower value we just skip this step and continue to the update.
1388 /* Select a new trial point.
1389 * If the derivatives at points a & c have different sign we interpolate to zero,
1390 * otherwise just do a bisection.
1392 if (gpa < 0 && gpc > 0)
1394 b = a + gpa*(a-c)/(gpc-gpa);
1401 /* safeguard if interpolation close to machine accuracy causes errors:
1402 * never go outside the interval
1404 if (b <= a || b >= c)
1409 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1411 /* Reload the old state */
1412 em_dd_partition_system(fplog, mdlog, -1, cr, top_global, inputrec,
1413 s_min, &top, mdAtoms, fr, vsite, constr,
1417 /* Take a trial step to this new point - new coords in s_b */
1418 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b,
1422 /* Calculate energy for the trial step */
1423 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1425 /* p does not change within a step, but since the domain decomposition
1426 * might change, we have to use cg_p of s_b here.
1428 const rvec *pb = s_b->s.cg_p.rvec_array();
1429 const rvec *sfb = s_b->f.rvec_array();
1431 for (int i = 0; i < mdatoms->homenr; i++)
1433 for (m = 0; m < DIM; m++)
1435 gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1438 /* Sum the gradient along the line across CPUs */
1441 gmx_sumd(1, &gpb, cr);
1446 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1447 s_a->epot, s_b->epot, s_c->epot, gpb);
1450 epot_repl = s_b->epot;
1452 /* Keep one of the intervals based on the value of the derivative at the new point */
1455 /* Replace c endpoint with b */
1456 swap_em_state(&s_b, &s_c);
1462 /* Replace a endpoint with b */
1463 swap_em_state(&s_b, &s_a);
1469 * Stop search as soon as we find a value smaller than the endpoints.
1470 * Never run more than 20 steps, no matter what.
1474 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1477 if (std::fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1480 /* OK. We couldn't find a significantly lower energy.
1481 * If beta==0 this was steepest descent, and then we give up.
1482 * If not, set beta=0 and restart with steepest descent before quitting.
1492 /* Reset memory before giving up */
1498 /* Select min energy state of A & C, put the best in B.
1500 if (s_c->epot < s_a->epot)
1504 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1505 s_c->epot, s_a->epot);
1507 swap_em_state(&s_b, &s_c);
1514 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1515 s_a->epot, s_c->epot);
1517 swap_em_state(&s_b, &s_a);
1526 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1529 swap_em_state(&s_b, &s_c);
1533 /* new search direction */
1534 /* beta = 0 means forget all memory and restart with steepest descents. */
1535 if (nstcg && ((step % nstcg) == 0))
1541 /* s_min->fnorm cannot be zero, because then we would have converged
1545 /* Polak-Ribiere update.
1546 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1548 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1550 /* Limit beta to prevent oscillations */
1551 if (fabs(beta) > 5.0)
1557 /* update positions */
1558 swap_em_state(&s_min, &s_b);
1561 /* Print it if necessary */
1564 if (mdrunOptions.verbose)
1566 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1567 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1568 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1569 s_min->fmax, s_min->a_fmax+1);
1572 /* Store the new (lower) energies */
1573 matrix nullBox = {};
1574 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(step),
1575 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1576 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1578 do_log = do_per_step(step, inputrec->nstlog);
1579 do_ene = do_per_step(step, inputrec->nstenergy);
1581 /* Prepare IMD energy record, if bIMD is TRUE. */
1582 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1586 print_ebin_header(fplog, step, step);
1588 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1589 do_log ? fplog : nullptr, step, step, eprNORMAL,
1590 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1593 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1594 if (MASTER(cr) && do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x.rvec_array(), inputrec, 0, wcycle))
1596 IMD_send_positions(inputrec->imd);
1599 /* Stop when the maximum force lies below tolerance.
1600 * If we have reached machine precision, converged is already set to true.
1602 converged = converged || (s_min->fmax < inputrec->em_tol);
1604 } /* End of the loop */
1606 /* IMD cleanup, if bIMD is TRUE. */
1607 IMD_finalize(inputrec->bIMD, inputrec->imd);
1611 step--; /* we never took that last step in this case */
1614 if (s_min->fmax > inputrec->em_tol)
1618 warn_step(fplog, inputrec->em_tol, s_min->fmax,
1619 step-1 == number_steps, FALSE);
1626 /* If we printed energy and/or logfile last step (which was the last step)
1627 * we don't have to do it again, but otherwise print the final values.
1631 /* Write final value to log since we didn't do anything the last step */
1632 print_ebin_header(fplog, step, step);
1634 if (!do_ene || !do_log)
1636 /* Write final energy file entries */
1637 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1638 !do_log ? fplog : nullptr, step, step, eprNORMAL,
1639 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1643 /* Print some stuff... */
1646 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1650 * For accurate normal mode calculation it is imperative that we
1651 * store the last conformation into the full precision binary trajectory.
1653 * However, we should only do it if we did NOT already write this step
1654 * above (which we did if do_x or do_f was true).
1656 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1657 * in the trajectory with the same step number.
1659 do_x = !do_per_step(step, inputrec->nstxout);
1660 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1662 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1663 top_global, inputrec, step,
1664 s_min, state_global, observablesHistory);
1669 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1670 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1671 s_min, sqrtNumAtoms);
1672 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1673 s_min, sqrtNumAtoms);
1675 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1678 finish_em(cr, outf, walltime_accounting, wcycle);
1680 /* To print the actual number of steps we needed somewhere */
1681 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1686 Integrator::do_lbfgs()
1688 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1691 gmx_enerdata_t *enerd;
1692 gmx_global_stat_t gstat;
1694 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1695 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1696 real *rho, *alpha, *p, *s, **dx, **dg;
1697 real a, b, c, maxdelta, delta;
1699 real dgdx, dgdg, sq, yr, beta;
1702 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1704 int start, end, number_steps;
1705 int i, k, m, n, gf, step;
1707 auto mdatoms = mdAtoms->mdatoms();
1709 GMX_LOG(mdlog.info).asParagraph().
1710 appendText("Note that activating L-BFGS energy minimization via the "
1711 "integrator .mdp option and the command gmx mdrun may "
1712 "be available in a different form in a future version of GROMACS, "
1713 "e.g. gmx minimize and an .mdp option.");
1717 gmx_fatal(FARGS, "L-BFGS minimization only supports a single rank");
1720 if (nullptr != constr)
1722 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).");
1725 n = 3*state_global->natoms;
1726 nmaxcorr = inputrec->nbfgscorr;
1731 snew(rho, nmaxcorr);
1732 snew(alpha, nmaxcorr);
1735 for (i = 0; i < nmaxcorr; i++)
1741 for (i = 0; i < nmaxcorr; i++)
1750 init_em(fplog, mdlog, LBFGS, cr, ms, inputrec, mdrunOptions,
1751 state_global, top_global, &ems, &top,
1752 nrnb, fr, &graph, mdAtoms, &gstat,
1753 vsite, constr, nullptr,
1755 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, inputrec, top_global, nullptr, wcycle);
1757 init_enerdata(top_global->groups.grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
1758 t_mdebin *mdebin = init_mdebin(mdoutf_get_fp_ene(outf), top_global, inputrec, nullptr);
1761 end = mdatoms->homenr;
1763 /* We need 4 working states */
1764 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1765 em_state_t *sa = &s0;
1766 em_state_t *sb = &s1;
1767 em_state_t *sc = &s2;
1768 em_state_t *last = &s3;
1769 /* Initialize by copying the state from ems (we could skip x and f here) */
1774 /* Print to log file */
1775 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1777 do_log = do_ene = do_x = do_f = TRUE;
1779 /* Max number of steps */
1780 number_steps = inputrec->nsteps;
1782 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1784 for (i = start; i < end; i++)
1786 if (mdatoms->cFREEZE)
1788 gf = mdatoms->cFREEZE[i];
1790 for (m = 0; m < DIM; m++)
1792 frozen[3*i+m] = (inputrec->opts.nFreeze[gf][m] != 0);
1797 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1801 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1806 construct_vsites(vsite, state_global->x.rvec_array(), 1, nullptr,
1807 top.idef.iparams, top.idef.il,
1808 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1811 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1812 /* do_force always puts the charge groups in the box and shifts again
1813 * We do not unshift, so molecules are always whole
1816 EnergyEvaluator energyEvaluator {
1817 fplog, mdlog, cr, ms,
1819 inputrec, nrnb, wcycle, gstat,
1820 vsite, constr, fcd, graph,
1821 mdAtoms, fr, ppForceWorkload, enerd
1823 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1827 /* Copy stuff to the energy bin for easy printing etc. */
1828 matrix nullBox = {};
1829 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(step),
1830 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
1831 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1833 print_ebin_header(fplog, step, step);
1834 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1835 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1838 /* Set the initial step.
1839 * since it will be multiplied by the non-normalized search direction
1840 * vector (force vector the first time), we scale it by the
1841 * norm of the force.
1846 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1847 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1848 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1849 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1850 fprintf(stderr, "\n");
1851 /* and copy to the log file too... */
1852 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1853 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1854 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1855 fprintf(fplog, "\n");
1858 // Point is an index to the memory of search directions, where 0 is the first one.
1861 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1862 real *fInit = static_cast<real *>(ems.f.rvec_array()[0]);
1863 for (i = 0; i < n; i++)
1867 dx[point][i] = fInit[i]; /* Initial search direction */
1875 // Stepsize will be modified during the search, and actually it is not critical
1876 // (the main efficiency in the algorithm comes from changing directions), but
1877 // we still need an initial value, so estimate it as the inverse of the norm
1878 // so we take small steps where the potential fluctuates a lot.
1879 stepsize = 1.0/ems.fnorm;
1881 /* Start the loop over BFGS steps.
1882 * Each successful step is counted, and we continue until
1883 * we either converge or reach the max number of steps.
1888 /* Set the gradient from the force */
1890 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1893 /* Write coordinates if necessary */
1894 do_x = do_per_step(step, inputrec->nstxout);
1895 do_f = do_per_step(step, inputrec->nstfout);
1900 mdof_flags |= MDOF_X;
1905 mdof_flags |= MDOF_F;
1910 mdof_flags |= MDOF_IMD;
1913 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1914 top_global, step, static_cast<real>(step), &ems.s, state_global, observablesHistory, ems.f);
1916 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1918 /* make s a pointer to current search direction - point=0 first time we get here */
1921 real *xx = static_cast<real *>(ems.s.x.rvec_array()[0]);
1922 real *ff = static_cast<real *>(ems.f.rvec_array()[0]);
1924 // calculate line gradient in position A
1925 for (gpa = 0, i = 0; i < n; i++)
1930 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1931 * relative change in coordinate is smaller than precision
1933 for (minstep = 0, i = 0; i < n; i++)
1943 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1945 if (stepsize < minstep)
1951 // Before taking any steps along the line, store the old position
1953 real *lastx = static_cast<real *>(last->s.x.data()[0]);
1954 real *lastf = static_cast<real *>(last->f.data()[0]);
1959 /* Take a step downhill.
1960 * In theory, we should find the actual minimum of the function in this
1961 * direction, somewhere along the line.
1962 * That is quite possible, but it turns out to take 5-10 function evaluations
1963 * for each line. However, we dont really need to find the exact minimum -
1964 * it is much better to start a new BFGS step in a modified direction as soon
1965 * as we are close to it. This will save a lot of energy evaluations.
1967 * In practice, we just try to take a single step.
1968 * If it worked (i.e. lowered the energy), we increase the stepsize but
1969 * continue straight to the next BFGS step without trying to find any minimum,
1970 * i.e. we change the search direction too. If the line was smooth, it is
1971 * likely we are in a smooth region, and then it makes sense to take longer
1972 * steps in the modified search direction too.
1974 * If it didn't work (higher energy), there must be a minimum somewhere between
1975 * the old position and the new one. Then we need to start by finding a lower
1976 * value before we change search direction. Since the energy was apparently
1977 * quite rough, we need to decrease the step size.
1979 * Due to the finite numerical accuracy, it turns out that it is a good idea
1980 * to accept a SMALL increase in energy, if the derivative is still downhill.
1981 * This leads to lower final energies in the tests I've done. / Erik
1984 // State "A" is the first position along the line.
1985 // reference position along line is initially zero
1988 // Check stepsize first. We do not allow displacements
1989 // larger than emstep.
1993 // Pick a new position C by adding stepsize to A.
1996 // Calculate what the largest change in any individual coordinate
1997 // would be (translation along line * gradient along line)
1999 for (i = 0; i < n; i++)
2002 if (delta > maxdelta)
2007 // If any displacement is larger than the stepsize limit, reduce the step
2008 if (maxdelta > inputrec->em_stepsize)
2013 while (maxdelta > inputrec->em_stepsize);
2015 // Take a trial step and move the coordinate array xc[] to position C
2016 real *xc = static_cast<real *>(sc->s.x.rvec_array()[0]);
2017 for (i = 0; i < n; i++)
2019 xc[i] = lastx[i] + c*s[i];
2023 // Calculate energy for the trial step in position C
2024 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
2026 // Calc line gradient in position C
2027 real *fc = static_cast<real *>(sc->f.rvec_array()[0]);
2028 for (gpc = 0, i = 0; i < n; i++)
2030 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2032 /* Sum the gradient along the line across CPUs */
2035 gmx_sumd(1, &gpc, cr);
2038 // This is the max amount of increase in energy we tolerate.
2039 // By allowing VERY small changes (close to numerical precision) we
2040 // frequently find even better (lower) final energies.
2041 tmp = std::sqrt(GMX_REAL_EPS)*fabs(sa->epot);
2043 // Accept the step if the energy is lower in the new position C (compared to A),
2044 // or if it is not significantly higher and the line derivative is still negative.
2045 foundlower = sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp));
2046 // If true, great, we found a better energy. We no longer try to alter the
2047 // stepsize, but simply accept this new better position. The we select a new
2048 // search direction instead, which will be much more efficient than continuing
2049 // to take smaller steps along a line. Set fnorm based on the new C position,
2050 // which will be used to update the stepsize to 1/fnorm further down.
2052 // If false, the energy is NOT lower in point C, i.e. it will be the same
2053 // or higher than in point A. In this case it is pointless to move to point C,
2054 // so we will have to do more iterations along the same line to find a smaller
2055 // value in the interval [A=0.0,C].
2056 // Here, A is still 0.0, but that will change when we do a search in the interval
2057 // [0.0,C] below. That search we will do by interpolation or bisection rather
2058 // than with the stepsize, so no need to modify it. For the next search direction
2059 // it will be reset to 1/fnorm anyway.
2063 // OK, if we didn't find a lower value we will have to locate one now - there must
2064 // be one in the interval [a,c].
2065 // The same thing is valid here, though: Don't spend dozens of iterations to find
2066 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2067 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2068 // I also have a safeguard for potentially really pathological functions so we never
2069 // take more than 20 steps before we give up.
2070 // If we already found a lower value we just skip this step and continue to the update.
2075 // Select a new trial point B in the interval [A,C].
2076 // If the derivatives at points a & c have different sign we interpolate to zero,
2077 // otherwise just do a bisection since there might be multiple minima/maxima
2078 // inside the interval.
2079 if (gpa < 0 && gpc > 0)
2081 b = a + gpa*(a-c)/(gpc-gpa);
2088 /* safeguard if interpolation close to machine accuracy causes errors:
2089 * never go outside the interval
2091 if (b <= a || b >= c)
2096 // Take a trial step to point B
2097 real *xb = static_cast<real *>(sb->s.x.rvec_array()[0]);
2098 for (i = 0; i < n; i++)
2100 xb[i] = lastx[i] + b*s[i];
2104 // Calculate energy for the trial step in point B
2105 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2108 // Calculate gradient in point B
2109 real *fb = static_cast<real *>(sb->f.rvec_array()[0]);
2110 for (gpb = 0, i = 0; i < n; i++)
2112 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2115 /* Sum the gradient along the line across CPUs */
2118 gmx_sumd(1, &gpb, cr);
2121 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2122 // at the new point B, and rename the endpoints of this new interval A and C.
2125 /* Replace c endpoint with b */
2127 /* swap states b and c */
2128 swap_em_state(&sb, &sc);
2132 /* Replace a endpoint with b */
2134 /* swap states a and b */
2135 swap_em_state(&sa, &sb);
2139 * Stop search as soon as we find a value smaller than the endpoints,
2140 * or if the tolerance is below machine precision.
2141 * Never run more than 20 steps, no matter what.
2145 while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2147 if (std::fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2149 /* OK. We couldn't find a significantly lower energy.
2150 * If ncorr==0 this was steepest descent, and then we give up.
2151 * If not, reset memory to restart as steepest descent before quitting.
2163 /* Search in gradient direction */
2164 for (i = 0; i < n; i++)
2166 dx[point][i] = ff[i];
2168 /* Reset stepsize */
2169 stepsize = 1.0/fnorm;
2174 /* Select min energy state of A & C, put the best in xx/ff/Epot
2176 if (sc->epot < sa->epot)
2198 /* Update the memory information, and calculate a new
2199 * approximation of the inverse hessian
2202 /* Have new data in Epot, xx, ff */
2203 if (ncorr < nmaxcorr)
2208 for (i = 0; i < n; i++)
2210 dg[point][i] = lastf[i]-ff[i];
2211 dx[point][i] *= step_taken;
2216 for (i = 0; i < n; i++)
2218 dgdg += dg[point][i]*dg[point][i];
2219 dgdx += dg[point][i]*dx[point][i];
2224 rho[point] = 1.0/dgdx;
2227 if (point >= nmaxcorr)
2233 for (i = 0; i < n; i++)
2240 /* Recursive update. First go back over the memory points */
2241 for (k = 0; k < ncorr; k++)
2250 for (i = 0; i < n; i++)
2252 sq += dx[cp][i]*p[i];
2255 alpha[cp] = rho[cp]*sq;
2257 for (i = 0; i < n; i++)
2259 p[i] -= alpha[cp]*dg[cp][i];
2263 for (i = 0; i < n; i++)
2268 /* And then go forward again */
2269 for (k = 0; k < ncorr; k++)
2272 for (i = 0; i < n; i++)
2274 yr += p[i]*dg[cp][i];
2278 beta = alpha[cp]-beta;
2280 for (i = 0; i < n; i++)
2282 p[i] += beta*dx[cp][i];
2292 for (i = 0; i < n; i++)
2296 dx[point][i] = p[i];
2304 /* Print it if necessary */
2307 if (mdrunOptions.verbose)
2309 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2310 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2311 step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2314 /* Store the new (lower) energies */
2315 matrix nullBox = {};
2316 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(step),
2317 mdatoms->tmass, enerd, nullptr, nullptr, nullptr, nullBox,
2318 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2319 do_log = do_per_step(step, inputrec->nstlog);
2320 do_ene = do_per_step(step, inputrec->nstenergy);
2323 print_ebin_header(fplog, step, step);
2325 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2326 do_log ? fplog : nullptr, step, step, eprNORMAL,
2327 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2330 /* Send x and E to IMD client, if bIMD is TRUE. */
2331 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x.rvec_array(), inputrec, 0, wcycle) && MASTER(cr))
2333 IMD_send_positions(inputrec->imd);
2336 // Reset stepsize in we are doing more iterations
2337 stepsize = 1.0/ems.fnorm;
2339 /* Stop when the maximum force lies below tolerance.
2340 * If we have reached machine precision, converged is already set to true.
2342 converged = converged || (ems.fmax < inputrec->em_tol);
2344 } /* End of the loop */
2346 /* IMD cleanup, if bIMD is TRUE. */
2347 IMD_finalize(inputrec->bIMD, inputrec->imd);
2351 step--; /* we never took that last step in this case */
2354 if (ems.fmax > inputrec->em_tol)
2358 warn_step(fplog, inputrec->em_tol, ems.fmax,
2359 step-1 == number_steps, FALSE);
2364 /* If we printed energy and/or logfile last step (which was the last step)
2365 * we don't have to do it again, but otherwise print the final values.
2367 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2369 print_ebin_header(fplog, step, step);
2371 if (!do_ene || !do_log) /* Write final energy file entries */
2373 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2374 !do_log ? fplog : nullptr, step, step, eprNORMAL,
2375 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2378 /* Print some stuff... */
2381 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2385 * For accurate normal mode calculation it is imperative that we
2386 * store the last conformation into the full precision binary trajectory.
2388 * However, we should only do it if we did NOT already write this step
2389 * above (which we did if do_x or do_f was true).
2391 do_x = !do_per_step(step, inputrec->nstxout);
2392 do_f = !do_per_step(step, inputrec->nstfout);
2393 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2394 top_global, inputrec, step,
2395 &ems, state_global, observablesHistory);
2399 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2400 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2401 number_steps, &ems, sqrtNumAtoms);
2402 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2403 number_steps, &ems, sqrtNumAtoms);
2405 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2408 finish_em(cr, outf, walltime_accounting, wcycle);
2410 /* To print the actual number of steps we needed somewhere */
2411 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2415 Integrator::do_steep()
2417 const char *SD = "Steepest Descents";
2419 gmx_enerdata_t *enerd;
2420 gmx_global_stat_t gstat;
2424 gmx_bool bDone, bAbort, do_x, do_f;
2429 int steps_accepted = 0;
2430 auto mdatoms = mdAtoms->mdatoms();
2432 GMX_LOG(mdlog.info).asParagraph().
2433 appendText("Note that activating steepest-descent energy minimization via the "
2434 "integrator .mdp option and the command gmx mdrun may "
2435 "be available in a different form in a future version of GROMACS, "
2436 "e.g. gmx minimize and an .mdp option.");
2438 /* Create 2 states on the stack and extract pointers that we will swap */
2439 em_state_t s0 {}, s1 {};
2440 em_state_t *s_min = &s0;
2441 em_state_t *s_try = &s1;
2443 /* Init em and store the local state in s_try */
2444 init_em(fplog, mdlog, SD, cr, ms, inputrec, mdrunOptions,
2445 state_global, top_global, s_try, &top,
2446 nrnb, fr, &graph, mdAtoms, &gstat,
2447 vsite, constr, nullptr,
2449 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, inputrec, top_global, nullptr, wcycle);
2451 init_enerdata(top_global->groups.grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
2452 t_mdebin *mdebin = init_mdebin(mdoutf_get_fp_ene(outf), top_global, inputrec, nullptr);
2454 /* Print to log file */
2455 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2457 /* Set variables for stepsize (in nm). This is the largest
2458 * step that we are going to make in any direction.
2460 ustep = inputrec->em_stepsize;
2463 /* Max number of steps */
2464 nsteps = inputrec->nsteps;
2468 /* Print to the screen */
2469 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2473 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2475 EnergyEvaluator energyEvaluator {
2476 fplog, mdlog, cr, ms,
2478 inputrec, nrnb, wcycle, gstat,
2479 vsite, constr, fcd, graph,
2480 mdAtoms, fr, ppForceWorkload, enerd
2483 /**** HERE STARTS THE LOOP ****
2484 * count is the counter for the number of steps
2485 * bDone will be TRUE when the minimization has converged
2486 * bAbort will be TRUE when nsteps steps have been performed or when
2487 * the stepsize becomes smaller than is reasonable for machine precision
2492 while (!bDone && !bAbort)
2494 bAbort = (nsteps >= 0) && (count == nsteps);
2496 /* set new coordinates, except for first step */
2497 bool validStep = true;
2501 do_em_step(cr, inputrec, mdatoms,
2502 s_min, stepsize, &s_min->f, s_try,
2508 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2512 // Signal constraint error during stepping with energy=inf
2513 s_try->epot = std::numeric_limits<real>::infinity();
2518 print_ebin_header(fplog, count, count);
2523 s_min->epot = s_try->epot;
2526 /* Print it if necessary */
2529 if (mdrunOptions.verbose)
2531 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2532 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2533 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2537 if ( (count == 0) || (s_try->epot < s_min->epot) )
2539 /* Store the new (lower) energies */
2540 matrix nullBox = {};
2541 upd_mdebin(mdebin, FALSE, FALSE, static_cast<double>(count),
2542 mdatoms->tmass, enerd, nullptr, nullptr, nullptr,
2543 nullBox, nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2545 /* Prepare IMD energy record, if bIMD is TRUE. */
2546 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2548 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2549 do_per_step(steps_accepted, inputrec->nstdisreout),
2550 do_per_step(steps_accepted, inputrec->nstorireout),
2551 fplog, count, count, eprNORMAL,
2552 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2557 /* Now if the new energy is smaller than the previous...
2558 * or if this is the first step!
2559 * or if we did random steps!
2562 if ( (count == 0) || (s_try->epot < s_min->epot) )
2566 /* Test whether the convergence criterion is met... */
2567 bDone = (s_try->fmax < inputrec->em_tol);
2569 /* Copy the arrays for force, positions and energy */
2570 /* The 'Min' array always holds the coords and forces of the minimal
2572 swap_em_state(&s_min, &s_try);
2578 /* Write to trn, if necessary */
2579 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2580 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2581 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2582 top_global, inputrec, count,
2583 s_min, state_global, observablesHistory);
2587 /* If energy is not smaller make the step smaller... */
2590 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2592 /* Reload the old state */
2593 em_dd_partition_system(fplog, mdlog, count, cr, top_global, inputrec,
2594 s_min, &top, mdAtoms, fr, vsite, constr,
2599 /* Determine new step */
2600 stepsize = ustep/s_min->fmax;
2602 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2604 if (count == nsteps || ustep < 1e-12)
2606 if (count == nsteps || ustep < 1e-6)
2611 warn_step(fplog, inputrec->em_tol, s_min->fmax,
2612 count == nsteps, constr != nullptr);
2617 /* Send IMD energies and positions, if bIMD is TRUE. */
2618 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box,
2619 MASTER(cr) ? state_global->x.rvec_array() : nullptr,
2620 inputrec, 0, wcycle) &&
2623 IMD_send_positions(inputrec->imd);
2627 } /* End of the loop */
2629 /* IMD cleanup, if bIMD is TRUE. */
2630 IMD_finalize(inputrec->bIMD, inputrec->imd);
2632 /* Print some data... */
2635 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2637 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout != 0, ftp2fn(efSTO, nfile, fnm),
2638 top_global, inputrec, count,
2639 s_min, state_global, observablesHistory);
2643 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2645 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2646 s_min, sqrtNumAtoms);
2647 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2648 s_min, sqrtNumAtoms);
2651 finish_em(cr, outf, walltime_accounting, wcycle);
2653 /* To print the actual number of steps we needed somewhere */
2654 inputrec->nsteps = count;
2656 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2662 const char *NM = "Normal Mode Analysis";
2665 gmx_enerdata_t *enerd;
2666 gmx_global_stat_t gstat;
2671 gmx_bool bSparse; /* use sparse matrix storage format */
2673 gmx_sparsematrix_t * sparse_matrix = nullptr;
2674 real * full_matrix = nullptr;
2676 /* added with respect to mdrun */
2678 real der_range = 10.0*std::sqrt(GMX_REAL_EPS);
2680 bool bIsMaster = MASTER(cr);
2681 auto mdatoms = mdAtoms->mdatoms();
2683 GMX_LOG(mdlog.info).asParagraph().
2684 appendText("Note that activating normal-mode analysis via the integrator "
2685 ".mdp option and the command gmx mdrun may "
2686 "be available in a different form in a future version of GROMACS, "
2687 "e.g. gmx normal-modes.");
2689 if (constr != nullptr)
2691 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2694 gmx_shellfc_t *shellfc;
2696 em_state_t state_work {};
2698 /* Init em and store the local state in state_minimum */
2699 init_em(fplog, mdlog, NM, cr, ms, inputrec, mdrunOptions,
2700 state_global, top_global, &state_work, &top,
2701 nrnb, fr, &graph, mdAtoms, &gstat,
2702 vsite, constr, &shellfc,
2704 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, inputrec, top_global, nullptr, wcycle);
2706 init_enerdata(top_global->groups.grps[egcENER].nr, inputrec->fepvals->n_lambda, enerd);
2708 std::vector<int> atom_index = get_atom_index(top_global);
2709 std::vector<gmx::RVec> fneg(atom_index.size(), {0, 0, 0});
2710 snew(dfdx, atom_index.size());
2716 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2717 " which MIGHT not be accurate enough for normal mode analysis.\n"
2718 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2719 " are fairly modest even if you recompile in double precision.\n\n");
2723 /* Check if we can/should use sparse storage format.
2725 * Sparse format is only useful when the Hessian itself is sparse, which it
2726 * will be when we use a cutoff.
2727 * For small systems (n<1000) it is easier to always use full matrix format, though.
2729 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2731 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2734 else if (atom_index.size() < 1000)
2736 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2742 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2746 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2747 sz = DIM*atom_index.size();
2749 fprintf(stderr, "Allocating Hessian memory...\n\n");
2753 sparse_matrix = gmx_sparsematrix_init(sz);
2754 sparse_matrix->compressed_symmetric = TRUE;
2758 snew(full_matrix, sz*sz);
2764 /* Write start time and temperature */
2765 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2767 /* fudge nr of steps to nr of atoms */
2768 inputrec->nsteps = atom_index.size()*2;
2772 fprintf(stderr, "starting normal mode calculation '%s'\n%" PRId64 " steps.\n\n",
2773 *(top_global->name), inputrec->nsteps);
2776 nnodes = cr->nnodes;
2778 /* Make evaluate_energy do a single node force calculation */
2780 EnergyEvaluator energyEvaluator {
2781 fplog, mdlog, cr, ms,
2783 inputrec, nrnb, wcycle, gstat,
2784 vsite, constr, fcd, graph,
2785 mdAtoms, fr, ppForceWorkload, enerd
2787 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2788 cr->nnodes = nnodes;
2790 /* if forces are not small, warn user */
2791 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2793 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2794 if (state_work.fmax > 1.0e-3)
2796 GMX_LOG(mdlog.warning).appendText(
2797 "The force is probably not small enough to "
2798 "ensure that you are at a minimum.\n"
2799 "Be aware that negative eigenvalues may occur\n"
2800 "when the resulting matrix is diagonalized.");
2803 /***********************************************************
2805 * Loop over all pairs in matrix
2807 * do_force called twice. Once with positive and
2808 * once with negative displacement
2810 ************************************************************/
2812 /* Steps are divided one by one over the nodes */
2814 auto state_work_x = makeArrayRef(state_work.s.x);
2815 auto state_work_f = makeArrayRef(state_work.f);
2816 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2818 size_t atom = atom_index[aid];
2819 for (size_t d = 0; d < DIM; d++)
2822 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2825 x_min = state_work_x[atom][d];
2827 for (unsigned int dx = 0; (dx < 2); dx++)
2831 state_work_x[atom][d] = x_min - der_range;
2835 state_work_x[atom][d] = x_min + der_range;
2838 /* Make evaluate_energy do a single node force calculation */
2842 /* Now is the time to relax the shells */
2843 relax_shell_flexcon(fplog,
2846 mdrunOptions.verbose,
2857 state_work.f.arrayRefWithPadding(),
2863 &top_global->groups,
2870 DdOpenBalanceRegionBeforeForceComputation::no,
2871 DdCloseBalanceRegionAfterForceComputation::no);
2877 energyEvaluator.run(&state_work, mu_tot, vir, pres, aid*2+dx, FALSE);
2880 cr->nnodes = nnodes;
2884 std::copy(state_work_f.begin(), state_work_f.begin()+atom_index.size(), fneg.begin());
2888 /* x is restored to original */
2889 state_work_x[atom][d] = x_min;
2891 for (size_t j = 0; j < atom_index.size(); j++)
2893 for (size_t k = 0; (k < DIM); k++)
2896 -(state_work_f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2903 #define mpi_type GMX_MPI_REAL
2904 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2905 cr->nodeid, cr->mpi_comm_mygroup);
2910 for (node = 0; (node < nnodes && aid+node < atom_index.size()); node++)
2916 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2917 cr->mpi_comm_mygroup, &stat);
2922 row = (aid + node)*DIM + d;
2924 for (size_t j = 0; j < atom_index.size(); j++)
2926 for (size_t k = 0; k < DIM; k++)
2932 if (col >= row && dfdx[j][k] != 0.0)
2934 gmx_sparsematrix_increment_value(sparse_matrix,
2935 row, col, dfdx[j][k]);
2940 full_matrix[row*sz+col] = dfdx[j][k];
2947 if (mdrunOptions.verbose && fplog)
2952 /* write progress */
2953 if (bIsMaster && mdrunOptions.verbose)
2955 fprintf(stderr, "\rFinished step %d out of %td",
2956 std::min<int>(atom+nnodes, atom_index.size()),
2964 fprintf(stderr, "\n\nWriting Hessian...\n");
2965 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2968 finish_em(cr, outf, walltime_accounting, wcycle);
2970 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);