2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 * \brief This file defines integrators for energy minimization
42 * \author Berk Hess <hess@kth.se>
43 * \author Erik Lindahl <erik@kth.se>
44 * \ingroup module_mdrun
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/collect.h"
59 #include "gromacs/domdec/dlbtiming.h"
60 #include "gromacs/domdec/domdec.h"
61 #include "gromacs/domdec/domdec_struct.h"
62 #include "gromacs/domdec/mdsetup.h"
63 #include "gromacs/domdec/partition.h"
64 #include "gromacs/ewald/pme_pp.h"
65 #include "gromacs/fileio/confio.h"
66 #include "gromacs/fileio/mtxio.h"
67 #include "gromacs/gmxlib/network.h"
68 #include "gromacs/gmxlib/nrnb.h"
69 #include "gromacs/imd/imd.h"
70 #include "gromacs/linearalgebra/sparsematrix.h"
71 #include "gromacs/listed_forces/manage_threading.h"
72 #include "gromacs/math/functions.h"
73 #include "gromacs/math/vec.h"
74 #include "gromacs/mdlib/constr.h"
75 #include "gromacs/mdlib/coupling.h"
76 #include "gromacs/mdlib/dispersioncorrection.h"
77 #include "gromacs/mdlib/ebin.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/energyoutput.h"
80 #include "gromacs/mdlib/force.h"
81 #include "gromacs/mdlib/force_flags.h"
82 #include "gromacs/mdlib/forcerec.h"
83 #include "gromacs/mdlib/gmx_omp_nthreads.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/stat.h"
87 #include "gromacs/mdlib/tgroup.h"
88 #include "gromacs/mdlib/trajectory_writing.h"
89 #include "gromacs/mdlib/update.h"
90 #include "gromacs/mdlib/vsite.h"
91 #include "gromacs/mdrunutility/handlerestart.h"
92 #include "gromacs/mdrunutility/printtime.h"
93 #include "gromacs/mdtypes/commrec.h"
94 #include "gromacs/mdtypes/forcerec.h"
95 #include "gromacs/mdtypes/inputrec.h"
96 #include "gromacs/mdtypes/interaction_const.h"
97 #include "gromacs/mdtypes/md_enums.h"
98 #include "gromacs/mdtypes/mdatom.h"
99 #include "gromacs/mdtypes/mdrunoptions.h"
100 #include "gromacs/mdtypes/state.h"
101 #include "gromacs/pbcutil/pbc.h"
102 #include "gromacs/timing/wallcycle.h"
103 #include "gromacs/timing/walltime_accounting.h"
104 #include "gromacs/topology/mtop_util.h"
105 #include "gromacs/topology/topology.h"
106 #include "gromacs/utility/cstringutil.h"
107 #include "gromacs/utility/exceptions.h"
108 #include "gromacs/utility/fatalerror.h"
109 #include "gromacs/utility/logger.h"
110 #include "gromacs/utility/smalloc.h"
112 #include "legacysimulator.h"
116 using gmx::MdrunScheduleWorkload
;
118 using gmx::VirtualSitesHandler
;
120 //! Utility structure for manipulating states during EM
121 typedef struct em_state
123 //! Copy of the global state
126 PaddedHostVector
<gmx::RVec
> f
;
129 //! Norm of the force
137 //! Print the EM starting conditions
138 static void print_em_start(FILE* fplog
,
140 gmx_walltime_accounting_t walltime_accounting
,
141 gmx_wallcycle_t wcycle
,
144 walltime_accounting_start_time(walltime_accounting
);
145 wallcycle_start(wcycle
, ewcRUN
);
146 print_start(fplog
, cr
, walltime_accounting
, name
);
149 //! Stop counting time for EM
150 static void em_time_end(gmx_walltime_accounting_t walltime_accounting
, gmx_wallcycle_t wcycle
)
152 wallcycle_stop(wcycle
, ewcRUN
);
154 walltime_accounting_end_time(walltime_accounting
);
157 //! Printing a log file and console header
158 static void sp_header(FILE* out
, const char* minimizer
, real ftol
, int nsteps
)
161 fprintf(out
, "%s:\n", minimizer
);
162 fprintf(out
, " Tolerance (Fmax) = %12.5e\n", ftol
);
163 fprintf(out
, " Number of steps = %12d\n", nsteps
);
166 //! Print warning message
167 static void warn_step(FILE* fp
, real ftol
, real fmax
, gmx_bool bLastStep
, gmx_bool bConstrain
)
169 constexpr bool realIsDouble
= GMX_DOUBLE
;
172 if (!std::isfinite(fmax
))
175 "\nEnergy minimization has stopped because the force "
176 "on at least one atom is not finite. This usually means "
177 "atoms are overlapping. Modify the input coordinates to "
178 "remove atom overlap or use soft-core potentials with "
179 "the free energy code to avoid infinite forces.\n%s",
180 !realIsDouble
? "You could also be lucky that switching to double precision "
181 "is sufficient to obtain finite forces.\n"
187 "\nEnergy minimization reached the maximum number "
188 "of steps before the forces reached the requested "
189 "precision Fmax < %g.\n",
195 "\nEnergy minimization has stopped, but the forces have "
196 "not converged to the requested precision Fmax < %g (which "
197 "may not be possible for your system). It stopped "
198 "because the algorithm tried to make a new step whose size "
199 "was too small, or there was no change in the energy since "
200 "last step. Either way, we regard the minimization as "
201 "converged to within the available machine precision, "
202 "given your starting configuration and EM parameters.\n%s%s",
204 !realIsDouble
? "\nDouble precision normally gives you higher accuracy, but "
205 "this is often not needed for preparing to run molecular "
208 bConstrain
? "You might need to increase your constraint accuracy, or turn\n"
209 "off constraints altogether (set constraints = none in mdp file)\n"
213 fputs(wrap_lines(buffer
, 78, 0, FALSE
), stderr
);
214 fputs(wrap_lines(buffer
, 78, 0, FALSE
), fp
);
217 //! Print message about convergence of the EM
218 static void print_converged(FILE* fp
,
224 const em_state_t
* ems
,
227 char buf
[STEPSTRSIZE
];
231 fprintf(fp
, "\n%s converged to Fmax < %g in %s steps\n", alg
, ftol
, gmx_step_str(count
, buf
));
233 else if (count
< nsteps
)
236 "\n%s converged to machine precision in %s steps,\n"
237 "but did not reach the requested Fmax < %g.\n",
238 alg
, gmx_step_str(count
, buf
), ftol
);
242 fprintf(fp
, "\n%s did not converge to Fmax < %g in %s steps.\n", alg
, ftol
,
243 gmx_step_str(count
, buf
));
247 fprintf(fp
, "Potential Energy = %21.14e\n", ems
->epot
);
248 fprintf(fp
, "Maximum force = %21.14e on atom %d\n", ems
->fmax
, ems
->a_fmax
+ 1);
249 fprintf(fp
, "Norm of force = %21.14e\n", ems
->fnorm
/ sqrtNumAtoms
);
251 fprintf(fp
, "Potential Energy = %14.7e\n", ems
->epot
);
252 fprintf(fp
, "Maximum force = %14.7e on atom %d\n", ems
->fmax
, ems
->a_fmax
+ 1);
253 fprintf(fp
, "Norm of force = %14.7e\n", ems
->fnorm
/ sqrtNumAtoms
);
257 //! Compute the norm and max of the force array in parallel
258 static void get_f_norm_max(const t_commrec
* cr
,
268 int la_max
, a_max
, start
, end
, i
, m
, gf
;
270 /* This routine finds the largest force and returns it.
271 * On parallel machines the global max is taken.
277 end
= mdatoms
->homenr
;
278 if (mdatoms
->cFREEZE
)
280 for (i
= start
; i
< end
; i
++)
282 gf
= mdatoms
->cFREEZE
[i
];
284 for (m
= 0; m
< DIM
; m
++)
286 if (!opts
->nFreeze
[gf
][m
])
288 fam
+= gmx::square(f
[i
][m
]);
301 for (i
= start
; i
< end
; i
++)
313 if (la_max
>= 0 && DOMAINDECOMP(cr
))
315 a_max
= cr
->dd
->globalAtomIndices
[la_max
];
323 snew(sum
, 2 * cr
->nnodes
+ 1);
324 sum
[2 * cr
->nodeid
] = fmax2
;
325 sum
[2 * cr
->nodeid
+ 1] = a_max
;
326 sum
[2 * cr
->nnodes
] = fnorm2
;
327 gmx_sumd(2 * cr
->nnodes
+ 1, sum
, cr
);
328 fnorm2
= sum
[2 * cr
->nnodes
];
329 /* Determine the global maximum */
330 for (i
= 0; i
< cr
->nnodes
; i
++)
332 if (sum
[2 * i
] > fmax2
)
335 a_max
= gmx::roundToInt(sum
[2 * i
+ 1]);
343 *fnorm
= sqrt(fnorm2
);
355 //! Compute the norm of the force
356 static void get_state_f_norm_max(const t_commrec
* cr
, t_grpopts
* opts
, t_mdatoms
* mdatoms
, em_state_t
* ems
)
358 get_f_norm_max(cr
, opts
, mdatoms
, ems
->f
.rvec_array(), &ems
->fnorm
, &ems
->fmax
, &ems
->a_fmax
);
361 //! Initialize the energy minimization
362 static void init_em(FILE* fplog
,
363 const gmx::MDLogger
& mdlog
,
367 gmx::ImdSession
* imdSession
,
369 t_state
* state_global
,
370 const gmx_mtop_t
* top_global
,
375 gmx::MDAtoms
* mdAtoms
,
376 gmx_global_stat_t
* gstat
,
377 VirtualSitesHandler
* vsite
,
378 gmx::Constraints
* constr
,
379 gmx_shellfc_t
** shellfc
)
385 fprintf(fplog
, "Initiating %s\n", title
);
390 state_global
->ngtc
= 0;
392 initialize_lambdas(fplog
, *ir
, MASTER(cr
), &(state_global
->fep_state
), state_global
->lambda
, nullptr);
396 GMX_ASSERT(shellfc
!= nullptr, "With NM we always support shells");
398 *shellfc
= init_shell_flexcon(stdout
, top_global
, constr
? constr
->numFlexibleConstraints() : 0,
399 ir
->nstcalcenergy
, DOMAINDECOMP(cr
));
403 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir
->eI
),
404 "This else currently only handles energy minimizers, consider if your algorithm "
405 "needs shell/flexible-constraint support");
407 /* With energy minimization, shells and flexible constraints are
408 * automatically minimized when treated like normal DOFS.
410 if (shellfc
!= nullptr)
416 if (DOMAINDECOMP(cr
))
418 dd_init_local_state(cr
->dd
, state_global
, &ems
->s
);
420 /* Distribute the charge groups over the nodes from the master node */
421 dd_partition_system(fplog
, mdlog
, ir
->init_step
, cr
, TRUE
, 1, state_global
, *top_global
, ir
,
422 imdSession
, pull_work
, &ems
->s
, &ems
->f
, mdAtoms
, top
, fr
, vsite
,
423 constr
, nrnb
, nullptr, FALSE
);
424 dd_store_state(cr
->dd
, &ems
->s
);
428 state_change_natoms(state_global
, state_global
->natoms
);
429 /* Just copy the state */
430 ems
->s
= *state_global
;
431 state_change_natoms(&ems
->s
, ems
->s
.natoms
);
433 mdAlgorithmsSetupAtomData(cr
, ir
, *top_global
, top
, fr
, &ems
->f
, mdAtoms
, constr
, vsite
,
434 shellfc
? *shellfc
: nullptr);
437 update_mdatoms(mdAtoms
->mdatoms(), ems
->s
.lambda
[efptMASS
]);
441 // TODO how should this cross-module support dependency be managed?
442 if (ir
->eConstrAlg
== econtSHAKE
&& gmx_mtop_ftype_count(top_global
, F_CONSTR
) > 0)
444 gmx_fatal(FARGS
, "Can not do energy minimization with %s, use %s\n",
445 econstr_names
[econtSHAKE
], econstr_names
[econtLINCS
]);
448 if (!ir
->bContinuation
)
450 /* Constrain the starting coordinates */
451 bool needsLogging
= true;
452 bool computeEnergy
= true;
453 bool computeVirial
= false;
455 constr
->apply(needsLogging
, computeEnergy
, -1, 0, 1.0, ems
->s
.x
.arrayRefWithPadding(),
456 ems
->s
.x
.arrayRefWithPadding(), ArrayRef
<RVec
>(), ems
->s
.box
,
457 ems
->s
.lambda
[efptFEP
], &dvdl_constr
, gmx::ArrayRefWithPadding
<RVec
>(),
458 computeVirial
, nullptr, gmx::ConstraintVariable::Positions
);
464 *gstat
= global_stat_init(ir
);
471 calc_shifts(ems
->s
.box
, fr
->shift_vec
);
474 //! Finalize the minimization
475 static void finish_em(const t_commrec
* cr
,
477 gmx_walltime_accounting_t walltime_accounting
,
478 gmx_wallcycle_t wcycle
)
480 if (!thisRankHasDuty(cr
, DUTY_PME
))
482 /* Tell the PME only node to finish */
483 gmx_pme_send_finish(cr
);
488 em_time_end(walltime_accounting
, wcycle
);
491 //! Swap two different EM states during minimization
492 static void swap_em_state(em_state_t
** ems1
, em_state_t
** ems2
)
501 //! Save the EM trajectory
502 static void write_em_traj(FILE* fplog
,
508 const gmx_mtop_t
* top_global
,
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
, top_global
->natoms
, step
,
533 static_cast<double>(step
), &state
->s
, state_global
,
534 observablesHistory
, state
->f
);
536 if (confout
!= nullptr)
538 if (DOMAINDECOMP(cr
))
540 /* If bX=true, x was collected to state_global in the call above */
543 auto globalXRef
= MASTER(cr
) ? state_global
->x
: gmx::ArrayRef
<gmx::RVec
>();
544 dd_collect_vec(cr
->dd
, &state
->s
, state
->s
.x
, globalXRef
);
549 /* Copy the local state pointer */
550 state_global
= &state
->s
;
555 if (ir
->pbcType
!= PbcType::No
&& !ir
->bPeriodicMols
&& DOMAINDECOMP(cr
))
557 /* Make molecules whole only for confout writing */
558 do_pbc_mtop(ir
->pbcType
, state
->s
.box
, top_global
, state_global
->x
.rvec_array());
561 write_sto_conf_mtop(confout
, *top_global
->name
, top_global
,
562 state_global
->x
.rvec_array(), nullptr, ir
->pbcType
, state
->s
.box
);
567 //! \brief Do one minimization step
569 // \returns true when the step succeeded, false when a constraint error occurred
570 static bool do_em_step(const t_commrec
* cr
,
575 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 /* OpenMP does not supported unsigned loop variables */
664 #pragma omp for schedule(static) nowait
665 for (gmx::index i
= 0; i
< gmx::ssize(s2
->cg_gl
); i
++)
667 s2
->cg_gl
[i
] = s1
->cg_gl
[i
];
672 if (DOMAINDECOMP(cr
))
674 s2
->ddp_count
= s1
->ddp_count
;
675 s2
->ddp_count_cg_gl
= s1
->ddp_count_cg_gl
;
681 validStep
= constr
->apply(
682 TRUE
, TRUE
, count
, 0, 1.0, s1
->x
.arrayRefWithPadding(), s2
->x
.arrayRefWithPadding(),
683 ArrayRef
<RVec
>(), s2
->box
, s2
->lambda
[efptBONDED
], &dvdl_constr
,
684 gmx::ArrayRefWithPadding
<RVec
>(), false, nullptr, gmx::ConstraintVariable::Positions
);
688 /* This global reduction will affect performance at high
689 * parallelization, but we can not really avoid it.
690 * But usually EM is not run at high parallelization.
692 int reductionBuffer
= static_cast<int>(!validStep
);
693 gmx_sumi(1, &reductionBuffer
, cr
);
694 validStep
= (reductionBuffer
== 0);
697 // We should move this check to the different minimizers
698 if (!validStep
&& ir
->eI
!= eiSteep
)
701 "The coordinates could not be constrained. Minimizer '%s' can not handle "
702 "constraint failures, use minimizer '%s' before using '%s'.",
703 EI(ir
->eI
), EI(eiSteep
), EI(ir
->eI
));
710 //! Prepare EM for using domain decomposition parallellization
711 static void em_dd_partition_system(FILE* fplog
,
712 const gmx::MDLogger
& mdlog
,
715 const gmx_mtop_t
* top_global
,
717 gmx::ImdSession
* imdSession
,
721 gmx::MDAtoms
* mdAtoms
,
723 VirtualSitesHandler
* vsite
,
724 gmx::Constraints
* constr
,
726 gmx_wallcycle_t wcycle
)
728 /* Repartition the domain decomposition */
729 dd_partition_system(fplog
, mdlog
, step
, cr
, FALSE
, 1, nullptr, *top_global
, ir
, imdSession
, pull_work
,
730 &ems
->s
, &ems
->f
, mdAtoms
, top
, fr
, vsite
, constr
, nrnb
, wcycle
, FALSE
);
731 dd_store_state(cr
->dd
, &ems
->s
);
737 /*! \brief Class to handle the work of setting and doing an energy evaluation.
739 * This class is a mere aggregate of parameters to pass to evaluate an
740 * energy, so that future changes to names and types of them consume
741 * less time when refactoring other code.
743 * Aggregate initialization is used, for which the chief risk is that
744 * if a member is added at the end and not all initializer lists are
745 * updated, then the member will be value initialized, which will
746 * typically mean initialization to zero.
748 * Use a braced initializer list to construct one of these. */
749 class EnergyEvaluator
752 /*! \brief Evaluates an energy on the state in \c ems.
754 * \todo In practice, the same objects mu_tot, vir, and pres
755 * are always passed to this function, so we would rather have
756 * them as data members. However, their C-array types are
757 * unsuited for aggregate initialization. When the types
758 * improve, the call signature of this method can be reduced.
760 void run(em_state_t
* ems
, rvec mu_tot
, tensor vir
, tensor pres
, int64_t count
, gmx_bool bFirst
);
761 //! Handles logging (deprecated).
764 const gmx::MDLogger
& mdlog
;
765 //! Handles communication.
767 //! Coordinates multi-simulations.
768 const gmx_multisim_t
* ms
;
769 //! Holds the simulation topology.
770 const gmx_mtop_t
* top_global
;
771 //! Holds the domain topology.
773 //! User input options.
774 t_inputrec
* inputrec
;
775 //! The Interactive Molecular Dynamics session.
776 gmx::ImdSession
* imdSession
;
777 //! The pull work object.
779 //! Manages flop accounting.
781 //! Manages wall cycle accounting.
782 gmx_wallcycle_t wcycle
;
783 //! Coordinates global reduction.
784 gmx_global_stat_t gstat
;
785 //! Handles virtual sites.
786 VirtualSitesHandler
* vsite
;
787 //! Handles constraints.
788 gmx::Constraints
* constr
;
789 //! Handles strange things.
791 //! Per-atom data for this domain.
792 gmx::MDAtoms
* mdAtoms
;
793 //! Handles how to calculate the forces.
795 //! Schedule of force-calculation work each step for this task.
796 MdrunScheduleWorkload
* runScheduleWork
;
797 //! Stores the computed energies.
798 gmx_enerdata_t
* enerd
;
801 void EnergyEvaluator::run(em_state_t
* ems
, rvec mu_tot
, tensor vir
, tensor pres
, int64_t count
, gmx_bool bFirst
)
805 tensor force_vir
, shake_vir
, ekin
;
809 /* Set the time to the initial time, the time does not change during EM */
810 t
= inputrec
->init_t
;
812 if (bFirst
|| (DOMAINDECOMP(cr
) && ems
->s
.ddp_count
< cr
->dd
->ddp_count
))
814 /* This is the first state or an old state used before the last ns */
820 if (inputrec
->nstlist
> 0)
828 vsite
->construct(ems
->s
.x
, 1, {}, ems
->s
.box
);
831 if (DOMAINDECOMP(cr
) && bNS
)
833 /* Repartition the domain decomposition */
834 em_dd_partition_system(fplog
, mdlog
, count
, cr
, top_global
, inputrec
, imdSession
, pull_work
,
835 ems
, top
, mdAtoms
, fr
, vsite
, constr
, nrnb
, wcycle
);
838 /* Calc force & energy on new trial position */
839 /* do_force always puts the charge groups in the box and shifts again
840 * We do not unshift, so molecules are always whole in congrad.c
842 do_force(fplog
, cr
, ms
, inputrec
, nullptr, nullptr, imdSession
, pull_work
, count
, nrnb
, wcycle
,
843 top
, ems
->s
.box
, ems
->s
.x
.arrayRefWithPadding(), &ems
->s
.hist
,
844 ems
->f
.arrayRefWithPadding(), force_vir
, mdAtoms
->mdatoms(), enerd
, fcd
, ems
->s
.lambda
,
845 fr
, runScheduleWork
, vsite
, mu_tot
, t
, nullptr,
846 GMX_FORCE_STATECHANGED
| GMX_FORCE_ALLFORCES
| GMX_FORCE_VIRIAL
| GMX_FORCE_ENERGY
847 | (bNS
? GMX_FORCE_NS
: 0),
848 DDBalanceRegionHandler(cr
));
850 /* Clear the unused shake virial and pressure */
851 clear_mat(shake_vir
);
854 /* Communicate stuff when parallel */
855 if (PAR(cr
) && inputrec
->eI
!= eiNM
)
857 wallcycle_start(wcycle
, ewcMoveE
);
859 global_stat(gstat
, cr
, enerd
, force_vir
, shake_vir
, inputrec
, nullptr, nullptr, nullptr, 1,
860 &terminate
, nullptr, FALSE
, CGLO_ENERGY
| CGLO_PRESSURE
| CGLO_CONSTRAINT
);
862 wallcycle_stop(wcycle
, ewcMoveE
);
865 if (fr
->dispersionCorrection
)
867 /* Calculate long range corrections to pressure and energy */
868 const DispersionCorrection::Correction correction
=
869 fr
->dispersionCorrection
->calculate(ems
->s
.box
, ems
->s
.lambda
[efptVDW
]);
871 enerd
->term
[F_DISPCORR
] = correction
.energy
;
872 enerd
->term
[F_EPOT
] += correction
.energy
;
873 enerd
->term
[F_PRES
] += correction
.pressure
;
874 enerd
->term
[F_DVDL
] += correction
.dvdl
;
878 enerd
->term
[F_DISPCORR
] = 0;
881 ems
->epot
= enerd
->term
[F_EPOT
];
885 /* Project out the constraint components of the force */
886 bool needsLogging
= false;
887 bool computeEnergy
= false;
888 bool computeVirial
= true;
890 auto f
= ems
->f
.arrayRefWithPadding();
891 constr
->apply(needsLogging
, computeEnergy
, count
, 0, 1.0, ems
->s
.x
.arrayRefWithPadding(), f
,
892 f
.unpaddedArrayRef(), ems
->s
.box
, ems
->s
.lambda
[efptBONDED
], &dvdl_constr
,
893 gmx::ArrayRefWithPadding
<RVec
>(), computeVirial
, shake_vir
,
894 gmx::ConstraintVariable::ForceDispl
);
895 enerd
->term
[F_DVDL_CONSTR
] += dvdl_constr
;
896 m_add(force_vir
, shake_vir
, vir
);
900 copy_mat(force_vir
, vir
);
904 enerd
->term
[F_PRES
] = calc_pres(fr
->pbcType
, inputrec
->nwall
, ems
->s
.box
, ekin
, vir
, pres
);
906 sum_dhdl(enerd
, ems
->s
.lambda
, *inputrec
->fepvals
);
908 if (EI_ENERGY_MINIMIZATION(inputrec
->eI
))
910 get_state_f_norm_max(cr
, &(inputrec
->opts
), mdAtoms
->mdatoms(), ems
);
916 //! Parallel utility summing energies and forces
917 static double reorder_partsum(const t_commrec
* cr
,
919 const gmx_mtop_t
* top_global
,
925 fprintf(debug
, "Doing reorder_partsum\n");
928 const rvec
* fm
= s_min
->f
.rvec_array();
929 const rvec
* fb
= s_b
->f
.rvec_array();
931 /* Collect fm in a global vector fmg.
932 * This conflicts with the spirit of domain decomposition,
933 * but to fully optimize this a much more complicated algorithm is required.
935 const int natoms
= top_global
->natoms
;
939 gmx::ArrayRef
<const int> indicesMin
= s_min
->s
.cg_gl
;
941 for (int a
: indicesMin
)
943 copy_rvec(fm
[i
], fmg
[a
]);
946 gmx_sum(top_global
->natoms
* 3, fmg
[0], cr
);
948 /* Now we will determine the part of the sum for the cgs in state s_b */
949 gmx::ArrayRef
<const int> indicesB
= s_b
->s
.cg_gl
;
954 gmx::ArrayRef
<const unsigned char> grpnrFREEZE
=
955 top_global
->groups
.groupNumbers
[SimulationAtomGroupType::Freeze
];
956 for (int a
: indicesB
)
958 if (!grpnrFREEZE
.empty())
962 for (int m
= 0; m
< DIM
; m
++)
964 if (!opts
->nFreeze
[gf
][m
])
966 partsum
+= (fb
[i
][m
] - fmg
[a
][m
]) * fb
[i
][m
];
977 //! Print some stuff, like beta, whatever that means.
978 static real
pr_beta(const t_commrec
* cr
,
981 const gmx_mtop_t
* top_global
,
987 /* This is just the classical Polak-Ribiere calculation of beta;
988 * it looks a bit complicated since we take freeze groups into account,
989 * and might have to sum it in parallel runs.
992 if (!DOMAINDECOMP(cr
)
993 || (s_min
->s
.ddp_count
== cr
->dd
->ddp_count
&& s_b
->s
.ddp_count
== cr
->dd
->ddp_count
))
995 const rvec
* fm
= s_min
->f
.rvec_array();
996 const rvec
* fb
= s_b
->f
.rvec_array();
999 /* This part of code can be incorrect with DD,
1000 * since the atom ordering in s_b and s_min might differ.
1002 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
1004 if (mdatoms
->cFREEZE
)
1006 gf
= mdatoms
->cFREEZE
[i
];
1008 for (int m
= 0; m
< DIM
; m
++)
1010 if (!opts
->nFreeze
[gf
][m
])
1012 sum
+= (fb
[i
][m
] - fm
[i
][m
]) * fb
[i
][m
];
1019 /* We need to reorder cgs while summing */
1020 sum
= reorder_partsum(cr
, opts
, top_global
, s_min
, s_b
);
1024 gmx_sumd(1, &sum
, cr
);
1027 return sum
/ gmx::square(s_min
->fnorm
);
1033 void LegacySimulator::do_cg()
1035 const char* CG
= "Polak-Ribiere Conjugate Gradients";
1037 gmx_localtop_t
top(top_global
->ffparams
);
1038 gmx_global_stat_t gstat
;
1039 double tmp
, minstep
;
1041 real a
, b
, c
, beta
= 0.0;
1044 gmx_bool converged
, foundlower
;
1045 rvec mu_tot
= { 0 };
1046 gmx_bool do_log
= FALSE
, do_ene
= FALSE
, do_x
, do_f
;
1048 int number_steps
, neval
= 0, nstcg
= inputrec
->nstcgsteep
;
1049 int m
, step
, nminstep
;
1050 auto mdatoms
= mdAtoms
->mdatoms();
1055 "Note that activating conjugate gradient energy minimization via the "
1056 "integrator .mdp option and the command gmx mdrun may "
1057 "be available in a different form in a future version of GROMACS, "
1058 "e.g. gmx minimize and an .mdp option.");
1064 // In CG, the state is extended with a search direction
1065 state_global
->flags
|= (1 << estCGP
);
1067 // Ensure the extra per-atom state array gets allocated
1068 state_change_natoms(state_global
, state_global
->natoms
);
1070 // Initialize the search direction to zero
1071 for (RVec
& cg_p
: state_global
->cg_p
)
1077 /* Create 4 states on the stack and extract pointers that we will swap */
1078 em_state_t s0
{}, s1
{}, s2
{}, s3
{};
1079 em_state_t
* s_min
= &s0
;
1080 em_state_t
* s_a
= &s1
;
1081 em_state_t
* s_b
= &s2
;
1082 em_state_t
* s_c
= &s3
;
1084 /* Init em and store the local state in s_min */
1085 init_em(fplog
, mdlog
, CG
, cr
, inputrec
, imdSession
, pull_work
, state_global
, top_global
, s_min
,
1086 &top
, nrnb
, fr
, mdAtoms
, &gstat
, vsite
, constr
, nullptr);
1087 const bool simulationsShareState
= false;
1088 gmx_mdoutf
* outf
= init_mdoutf(fplog
, nfile
, fnm
, mdrunOptions
, cr
, outputProvider
,
1089 mdModulesNotifier
, inputrec
, top_global
, nullptr, wcycle
,
1090 StartingBehavior::NewSimulation
, simulationsShareState
, ms
);
1091 gmx::EnergyOutput
energyOutput(mdoutf_get_fp_ene(outf
), top_global
, inputrec
, pull_work
, nullptr,
1092 false, StartingBehavior::NewSimulation
, mdModulesNotifier
);
1094 /* Print to log file */
1095 print_em_start(fplog
, cr
, walltime_accounting
, wcycle
, CG
);
1097 /* Max number of steps */
1098 number_steps
= inputrec
->nsteps
;
1102 sp_header(stderr
, CG
, inputrec
->em_tol
, number_steps
);
1106 sp_header(fplog
, CG
, inputrec
->em_tol
, number_steps
);
1109 EnergyEvaluator energyEvaluator
{
1110 fplog
, mdlog
, cr
, ms
, top_global
, &top
, inputrec
, imdSession
, pull_work
,
1111 nrnb
, wcycle
, gstat
, vsite
, constr
, fcd
, mdAtoms
, fr
, runScheduleWork
,
1114 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1115 /* do_force always puts the charge groups in the box and shifts again
1116 * We do not unshift, so molecules are always whole in congrad.c
1118 energyEvaluator
.run(s_min
, mu_tot
, vir
, pres
, -1, TRUE
);
1122 /* Copy stuff to the energy bin for easy printing etc. */
1123 matrix nullBox
= {};
1124 energyOutput
.addDataAtEnergyStep(false, false, static_cast<double>(step
), mdatoms
->tmass
,
1125 enerd
, nullptr, nullptr, nullptr, nullBox
, nullptr,
1126 nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
1128 EnergyOutput::printHeader(fplog
, step
, step
);
1129 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), TRUE
, FALSE
, FALSE
, fplog
, step
,
1130 step
, fcd
, nullptr);
1133 /* Estimate/guess the initial stepsize */
1134 stepsize
= inputrec
->em_stepsize
/ s_min
->fnorm
;
1138 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
1139 fprintf(stderr
, " F-max = %12.5e on atom %d\n", s_min
->fmax
, s_min
->a_fmax
+ 1);
1140 fprintf(stderr
, " F-Norm = %12.5e\n", s_min
->fnorm
/ sqrtNumAtoms
);
1141 fprintf(stderr
, "\n");
1142 /* and copy to the log file too... */
1143 fprintf(fplog
, " F-max = %12.5e on atom %d\n", s_min
->fmax
, s_min
->a_fmax
+ 1);
1144 fprintf(fplog
, " F-Norm = %12.5e\n", s_min
->fnorm
/ sqrtNumAtoms
);
1145 fprintf(fplog
, "\n");
1147 /* Start the loop over CG steps.
1148 * Each successful step is counted, and we continue until
1149 * we either converge or reach the max number of steps.
1152 for (step
= 0; (number_steps
< 0 || step
<= number_steps
) && !converged
; step
++)
1155 /* start taking steps in a new direction
1156 * First time we enter the routine, beta=0, and the direction is
1157 * simply the negative gradient.
1160 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1161 rvec
* pm
= s_min
->s
.cg_p
.rvec_array();
1162 const rvec
* sfm
= s_min
->f
.rvec_array();
1165 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
1167 if (mdatoms
->cFREEZE
)
1169 gf
= mdatoms
->cFREEZE
[i
];
1171 for (m
= 0; m
< DIM
; m
++)
1173 if (!inputrec
->opts
.nFreeze
[gf
][m
])
1175 pm
[i
][m
] = sfm
[i
][m
] + beta
* pm
[i
][m
];
1176 gpa
-= pm
[i
][m
] * sfm
[i
][m
];
1177 /* f is negative gradient, thus the sign */
1186 /* Sum the gradient along the line across CPUs */
1189 gmx_sumd(1, &gpa
, cr
);
1192 /* Calculate the norm of the search vector */
1193 get_f_norm_max(cr
, &(inputrec
->opts
), mdatoms
, pm
, &pnorm
, nullptr, nullptr);
1195 /* Just in case stepsize reaches zero due to numerical precision... */
1198 stepsize
= inputrec
->em_stepsize
/ pnorm
;
1202 * Double check the value of the derivative in the search direction.
1203 * If it is positive it must be due to the old information in the
1204 * CG formula, so just remove that and start over with beta=0.
1205 * This corresponds to a steepest descent step.
1210 step
--; /* Don't count this step since we are restarting */
1211 continue; /* Go back to the beginning of the big for-loop */
1214 /* Calculate minimum allowed stepsize, before the average (norm)
1215 * relative change in coordinate is smaller than precision
1218 auto s_min_x
= makeArrayRef(s_min
->s
.x
);
1219 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
1221 for (m
= 0; m
< DIM
; m
++)
1223 tmp
= fabs(s_min_x
[i
][m
]);
1228 tmp
= pm
[i
][m
] / tmp
;
1229 minstep
+= tmp
* tmp
;
1232 /* Add up from all CPUs */
1235 gmx_sumd(1, &minstep
, cr
);
1238 minstep
= GMX_REAL_EPS
/ sqrt(minstep
/ (3 * top_global
->natoms
));
1240 if (stepsize
< minstep
)
1246 /* Write coordinates if necessary */
1247 do_x
= do_per_step(step
, inputrec
->nstxout
);
1248 do_f
= do_per_step(step
, inputrec
->nstfout
);
1250 write_em_traj(fplog
, cr
, outf
, do_x
, do_f
, nullptr, top_global
, inputrec
, step
, s_min
,
1251 state_global
, observablesHistory
);
1253 /* Take a step downhill.
1254 * In theory, we should minimize the function along this direction.
1255 * That is quite possible, but it turns out to take 5-10 function evaluations
1256 * for each line. However, we dont really need to find the exact minimum -
1257 * it is much better to start a new CG step in a modified direction as soon
1258 * as we are close to it. This will save a lot of energy evaluations.
1260 * In practice, we just try to take a single step.
1261 * If it worked (i.e. lowered the energy), we increase the stepsize but
1262 * the continue straight to the next CG step without trying to find any minimum.
1263 * If it didn't work (higher energy), there must be a minimum somewhere between
1264 * the old position and the new one.
1266 * Due to the finite numerical accuracy, it turns out that it is a good idea
1267 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1268 * This leads to lower final energies in the tests I've done. / Erik
1270 s_a
->epot
= s_min
->epot
;
1272 c
= a
+ stepsize
; /* reference position along line is zero */
1274 if (DOMAINDECOMP(cr
) && s_min
->s
.ddp_count
< cr
->dd
->ddp_count
)
1276 em_dd_partition_system(fplog
, mdlog
, step
, cr
, top_global
, inputrec
, imdSession
,
1277 pull_work
, s_min
, &top
, mdAtoms
, fr
, vsite
, constr
, nrnb
, wcycle
);
1280 /* Take a trial step (new coords in s_c) */
1281 do_em_step(cr
, inputrec
, mdatoms
, s_min
, c
, &s_min
->s
.cg_p
, s_c
, constr
, -1);
1284 /* Calculate energy for the trial step */
1285 energyEvaluator
.run(s_c
, mu_tot
, vir
, pres
, -1, FALSE
);
1287 /* Calc derivative along line */
1288 const rvec
* pc
= s_c
->s
.cg_p
.rvec_array();
1289 const rvec
* sfc
= s_c
->f
.rvec_array();
1291 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
1293 for (m
= 0; m
< DIM
; m
++)
1295 gpc
-= pc
[i
][m
] * sfc
[i
][m
]; /* f is negative gradient, thus the sign */
1298 /* Sum the gradient along the line across CPUs */
1301 gmx_sumd(1, &gpc
, cr
);
1304 /* This is the max amount of increase in energy we tolerate */
1305 tmp
= std::sqrt(GMX_REAL_EPS
) * fabs(s_a
->epot
);
1307 /* Accept the step if the energy is lower, or if it is not significantly higher
1308 * and the line derivative is still negative.
1310 if (s_c
->epot
< s_a
->epot
|| (gpc
< 0 && s_c
->epot
< (s_a
->epot
+ tmp
)))
1313 /* Great, we found a better energy. Increase step for next iteration
1314 * if we are still going down, decrease it otherwise
1318 stepsize
*= 1.618034; /* The golden section */
1322 stepsize
*= 0.618034; /* 1/golden section */
1327 /* New energy is the same or higher. We will have to do some work
1328 * to find a smaller value in the interval. Take smaller step next time!
1331 stepsize
*= 0.618034;
1335 /* OK, if we didn't find a lower value we will have to locate one now - there must
1336 * be one in the interval [a=0,c].
1337 * The same thing is valid here, though: Don't spend dozens of iterations to find
1338 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1339 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1341 * I also have a safeguard for potentially really pathological functions so we never
1342 * take more than 20 steps before we give up ...
1344 * If we already found a lower value we just skip this step and continue to the update.
1353 /* Select a new trial point.
1354 * If the derivatives at points a & c have different sign we interpolate to zero,
1355 * otherwise just do a bisection.
1357 if (gpa
< 0 && gpc
> 0)
1359 b
= a
+ gpa
* (a
- c
) / (gpc
- gpa
);
1366 /* safeguard if interpolation close to machine accuracy causes errors:
1367 * never go outside the interval
1369 if (b
<= a
|| b
>= c
)
1374 if (DOMAINDECOMP(cr
) && s_min
->s
.ddp_count
!= cr
->dd
->ddp_count
)
1376 /* Reload the old state */
1377 em_dd_partition_system(fplog
, mdlog
, -1, cr
, top_global
, inputrec
, imdSession
, pull_work
,
1378 s_min
, &top
, mdAtoms
, fr
, vsite
, constr
, nrnb
, wcycle
);
1381 /* Take a trial step to this new point - new coords in s_b */
1382 do_em_step(cr
, inputrec
, mdatoms
, s_min
, b
, &s_min
->s
.cg_p
, s_b
, constr
, -1);
1385 /* Calculate energy for the trial step */
1386 energyEvaluator
.run(s_b
, mu_tot
, vir
, pres
, -1, FALSE
);
1388 /* p does not change within a step, but since the domain decomposition
1389 * might change, we have to use cg_p of s_b here.
1391 const rvec
* pb
= s_b
->s
.cg_p
.rvec_array();
1392 const rvec
* sfb
= s_b
->f
.rvec_array();
1394 for (int i
= 0; i
< mdatoms
->homenr
; i
++)
1396 for (m
= 0; m
< DIM
; m
++)
1398 gpb
-= pb
[i
][m
] * sfb
[i
][m
]; /* f is negative gradient, thus the sign */
1401 /* Sum the gradient along the line across CPUs */
1404 gmx_sumd(1, &gpb
, cr
);
1409 fprintf(debug
, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n", s_a
->epot
, s_b
->epot
,
1413 epot_repl
= s_b
->epot
;
1415 /* Keep one of the intervals based on the value of the derivative at the new point */
1418 /* Replace c endpoint with b */
1419 swap_em_state(&s_b
, &s_c
);
1425 /* Replace a endpoint with b */
1426 swap_em_state(&s_b
, &s_a
);
1432 * Stop search as soon as we find a value smaller than the endpoints.
1433 * Never run more than 20 steps, no matter what.
1436 } while ((epot_repl
> s_a
->epot
|| epot_repl
> s_c
->epot
) && (nminstep
< 20));
1438 if (std::fabs(epot_repl
- s_min
->epot
) < fabs(s_min
->epot
) * GMX_REAL_EPS
|| nminstep
>= 20)
1440 /* OK. We couldn't find a significantly lower energy.
1441 * If beta==0 this was steepest descent, and then we give up.
1442 * If not, set beta=0 and restart with steepest descent before quitting.
1452 /* Reset memory before giving up */
1458 /* Select min energy state of A & C, put the best in B.
1460 if (s_c
->epot
< s_a
->epot
)
1464 fprintf(debug
, "CGE: C (%f) is lower than A (%f), moving C to B\n", s_c
->epot
,
1467 swap_em_state(&s_b
, &s_c
);
1474 fprintf(debug
, "CGE: A (%f) is lower than C (%f), moving A to B\n", s_a
->epot
,
1477 swap_em_state(&s_b
, &s_a
);
1485 fprintf(debug
, "CGE: Found a lower energy %f, moving C to B\n", s_c
->epot
);
1487 swap_em_state(&s_b
, &s_c
);
1491 /* new search direction */
1492 /* beta = 0 means forget all memory and restart with steepest descents. */
1493 if (nstcg
&& ((step
% nstcg
) == 0))
1499 /* s_min->fnorm cannot be zero, because then we would have converged
1503 /* Polak-Ribiere update.
1504 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1506 beta
= pr_beta(cr
, &inputrec
->opts
, mdatoms
, top_global
, s_min
, s_b
);
1508 /* Limit beta to prevent oscillations */
1509 if (fabs(beta
) > 5.0)
1515 /* update positions */
1516 swap_em_state(&s_min
, &s_b
);
1519 /* Print it if necessary */
1522 if (mdrunOptions
.verbose
)
1524 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
1525 fprintf(stderr
, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step
,
1526 s_min
->epot
, s_min
->fnorm
/ sqrtNumAtoms
, s_min
->fmax
, s_min
->a_fmax
+ 1);
1529 /* Store the new (lower) energies */
1530 matrix nullBox
= {};
1531 energyOutput
.addDataAtEnergyStep(false, false, static_cast<double>(step
), mdatoms
->tmass
,
1532 enerd
, nullptr, nullptr, nullptr, nullBox
, nullptr,
1533 nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
1535 do_log
= do_per_step(step
, inputrec
->nstlog
);
1536 do_ene
= do_per_step(step
, inputrec
->nstenergy
);
1538 imdSession
->fillEnergyRecord(step
, TRUE
);
1542 EnergyOutput::printHeader(fplog
, step
, step
);
1544 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), do_ene
, FALSE
, FALSE
,
1545 do_log
? fplog
: nullptr, step
, step
, fcd
, nullptr);
1548 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1549 if (MASTER(cr
) && imdSession
->run(step
, TRUE
, state_global
->box
, state_global
->x
.rvec_array(), 0))
1551 imdSession
->sendPositionsAndEnergies();
1554 /* Stop when the maximum force lies below tolerance.
1555 * If we have reached machine precision, converged is already set to true.
1557 converged
= converged
|| (s_min
->fmax
< inputrec
->em_tol
);
1559 } /* End of the loop */
1563 step
--; /* we never took that last step in this case */
1565 if (s_min
->fmax
> inputrec
->em_tol
)
1569 warn_step(fplog
, inputrec
->em_tol
, s_min
->fmax
, step
- 1 == number_steps
, FALSE
);
1576 /* If we printed energy and/or logfile last step (which was the last step)
1577 * we don't have to do it again, but otherwise print the final values.
1581 /* Write final value to log since we didn't do anything the last step */
1582 EnergyOutput::printHeader(fplog
, step
, step
);
1584 if (!do_ene
|| !do_log
)
1586 /* Write final energy file entries */
1587 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), !do_ene
, FALSE
, FALSE
,
1588 !do_log
? fplog
: nullptr, step
, step
, fcd
, nullptr);
1592 /* Print some stuff... */
1595 fprintf(stderr
, "\nwriting lowest energy coordinates.\n");
1599 * For accurate normal mode calculation it is imperative that we
1600 * store the last conformation into the full precision binary trajectory.
1602 * However, we should only do it if we did NOT already write this step
1603 * above (which we did if do_x or do_f was true).
1605 /* Note that with 0 < nstfout != nstxout we can end up with two frames
1606 * in the trajectory with the same step number.
1608 do_x
= !do_per_step(step
, inputrec
->nstxout
);
1609 do_f
= (inputrec
->nstfout
> 0 && !do_per_step(step
, inputrec
->nstfout
));
1611 write_em_traj(fplog
, cr
, outf
, do_x
, do_f
, ftp2fn(efSTO
, nfile
, fnm
), top_global
, inputrec
,
1612 step
, s_min
, state_global
, observablesHistory
);
1617 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
1618 print_converged(stderr
, CG
, inputrec
->em_tol
, step
, converged
, number_steps
, s_min
, sqrtNumAtoms
);
1619 print_converged(fplog
, CG
, inputrec
->em_tol
, step
, converged
, number_steps
, s_min
, sqrtNumAtoms
);
1621 fprintf(fplog
, "\nPerformed %d energy evaluations in total.\n", neval
);
1624 finish_em(cr
, outf
, walltime_accounting
, wcycle
);
1626 /* To print the actual number of steps we needed somewhere */
1627 walltime_accounting_set_nsteps_done(walltime_accounting
, step
);
1631 void LegacySimulator::do_lbfgs()
1633 static const char* LBFGS
= "Low-Memory BFGS Minimizer";
1635 gmx_localtop_t
top(top_global
->ffparams
);
1636 gmx_global_stat_t gstat
;
1637 int ncorr
, nmaxcorr
, point
, cp
, neval
, nminstep
;
1638 double stepsize
, step_taken
, gpa
, gpb
, gpc
, tmp
, minstep
;
1639 real
* rho
, *alpha
, *p
, *s
, **dx
, **dg
;
1640 real a
, b
, c
, maxdelta
, delta
;
1642 real dgdx
, dgdg
, sq
, yr
, beta
;
1644 rvec mu_tot
= { 0 };
1645 gmx_bool do_log
, do_ene
, do_x
, do_f
, foundlower
, *frozen
;
1647 int start
, end
, number_steps
;
1648 int i
, k
, m
, n
, gf
, step
;
1650 auto mdatoms
= mdAtoms
->mdatoms();
1655 "Note that activating L-BFGS energy minimization via the "
1656 "integrator .mdp option and the command gmx mdrun may "
1657 "be available in a different form in a future version of GROMACS, "
1658 "e.g. gmx minimize and an .mdp option.");
1662 gmx_fatal(FARGS
, "L-BFGS minimization only supports a single rank");
1665 if (nullptr != constr
)
1669 "The combination of constraints and L-BFGS minimization is not implemented. Either "
1670 "do not use constraints, or use another minimizer (e.g. steepest descent).");
1673 n
= 3 * state_global
->natoms
;
1674 nmaxcorr
= inputrec
->nbfgscorr
;
1679 snew(rho
, nmaxcorr
);
1680 snew(alpha
, nmaxcorr
);
1683 for (i
= 0; i
< nmaxcorr
; i
++)
1689 for (i
= 0; i
< nmaxcorr
; i
++)
1698 init_em(fplog
, mdlog
, LBFGS
, cr
, inputrec
, imdSession
, pull_work
, state_global
, top_global
,
1699 &ems
, &top
, nrnb
, fr
, mdAtoms
, &gstat
, vsite
, constr
, nullptr);
1700 const bool simulationsShareState
= false;
1701 gmx_mdoutf
* outf
= init_mdoutf(fplog
, nfile
, fnm
, mdrunOptions
, cr
, outputProvider
,
1702 mdModulesNotifier
, inputrec
, top_global
, nullptr, wcycle
,
1703 StartingBehavior::NewSimulation
, simulationsShareState
, ms
);
1704 gmx::EnergyOutput
energyOutput(mdoutf_get_fp_ene(outf
), top_global
, inputrec
, pull_work
, nullptr,
1705 false, StartingBehavior::NewSimulation
, mdModulesNotifier
);
1708 end
= mdatoms
->homenr
;
1710 /* We need 4 working states */
1711 em_state_t s0
{}, s1
{}, s2
{}, s3
{};
1712 em_state_t
* sa
= &s0
;
1713 em_state_t
* sb
= &s1
;
1714 em_state_t
* sc
= &s2
;
1715 em_state_t
* last
= &s3
;
1716 /* Initialize by copying the state from ems (we could skip x and f here) */
1721 /* Print to log file */
1722 print_em_start(fplog
, cr
, walltime_accounting
, wcycle
, LBFGS
);
1724 do_log
= do_ene
= do_x
= do_f
= TRUE
;
1726 /* Max number of steps */
1727 number_steps
= inputrec
->nsteps
;
1729 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1731 for (i
= start
; i
< end
; i
++)
1733 if (mdatoms
->cFREEZE
)
1735 gf
= mdatoms
->cFREEZE
[i
];
1737 for (m
= 0; m
< DIM
; m
++)
1739 frozen
[3 * i
+ m
] = (inputrec
->opts
.nFreeze
[gf
][m
] != 0);
1744 sp_header(stderr
, LBFGS
, inputrec
->em_tol
, number_steps
);
1748 sp_header(fplog
, LBFGS
, inputrec
->em_tol
, number_steps
);
1753 vsite
->construct(state_global
->x
, 1, {}, state_global
->box
);
1756 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1757 /* do_force always puts the charge groups in the box and shifts again
1758 * We do not unshift, so molecules are always whole
1761 EnergyEvaluator energyEvaluator
{
1762 fplog
, mdlog
, cr
, ms
, top_global
, &top
, inputrec
, imdSession
, pull_work
,
1763 nrnb
, wcycle
, gstat
, vsite
, constr
, fcd
, mdAtoms
, fr
, runScheduleWork
,
1766 energyEvaluator
.run(&ems
, mu_tot
, vir
, pres
, -1, TRUE
);
1770 /* Copy stuff to the energy bin for easy printing etc. */
1771 matrix nullBox
= {};
1772 energyOutput
.addDataAtEnergyStep(false, false, static_cast<double>(step
), mdatoms
->tmass
,
1773 enerd
, nullptr, nullptr, nullptr, nullBox
, nullptr,
1774 nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
1776 EnergyOutput::printHeader(fplog
, step
, step
);
1777 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), TRUE
, FALSE
, FALSE
, fplog
, step
,
1778 step
, fcd
, nullptr);
1781 /* Set the initial step.
1782 * since it will be multiplied by the non-normalized search direction
1783 * vector (force vector the first time), we scale it by the
1784 * norm of the force.
1789 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
1790 fprintf(stderr
, "Using %d BFGS correction steps.\n\n", nmaxcorr
);
1791 fprintf(stderr
, " F-max = %12.5e on atom %d\n", ems
.fmax
, ems
.a_fmax
+ 1);
1792 fprintf(stderr
, " F-Norm = %12.5e\n", ems
.fnorm
/ sqrtNumAtoms
);
1793 fprintf(stderr
, "\n");
1794 /* and copy to the log file too... */
1795 fprintf(fplog
, "Using %d BFGS correction steps.\n\n", nmaxcorr
);
1796 fprintf(fplog
, " F-max = %12.5e on atom %d\n", ems
.fmax
, ems
.a_fmax
+ 1);
1797 fprintf(fplog
, " F-Norm = %12.5e\n", ems
.fnorm
/ sqrtNumAtoms
);
1798 fprintf(fplog
, "\n");
1801 // Point is an index to the memory of search directions, where 0 is the first one.
1804 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1805 real
* fInit
= static_cast<real
*>(ems
.f
.rvec_array()[0]);
1806 for (i
= 0; i
< n
; i
++)
1810 dx
[point
][i
] = fInit
[i
]; /* Initial search direction */
1818 // Stepsize will be modified during the search, and actually it is not critical
1819 // (the main efficiency in the algorithm comes from changing directions), but
1820 // we still need an initial value, so estimate it as the inverse of the norm
1821 // so we take small steps where the potential fluctuates a lot.
1822 stepsize
= 1.0 / ems
.fnorm
;
1824 /* Start the loop over BFGS steps.
1825 * Each successful step is counted, and we continue until
1826 * we either converge or reach the max number of steps.
1831 /* Set the gradient from the force */
1833 for (step
= 0; (number_steps
< 0 || step
<= number_steps
) && !converged
; step
++)
1836 /* Write coordinates if necessary */
1837 do_x
= do_per_step(step
, inputrec
->nstxout
);
1838 do_f
= do_per_step(step
, inputrec
->nstfout
);
1843 mdof_flags
|= MDOF_X
;
1848 mdof_flags
|= MDOF_F
;
1853 mdof_flags
|= MDOF_IMD
;
1856 mdoutf_write_to_trajectory_files(fplog
, cr
, outf
, mdof_flags
, top_global
->natoms
, step
,
1857 static_cast<real
>(step
), &ems
.s
, state_global
,
1858 observablesHistory
, ems
.f
);
1860 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1862 /* make s a pointer to current search direction - point=0 first time we get here */
1865 real
* xx
= static_cast<real
*>(ems
.s
.x
.rvec_array()[0]);
1866 real
* ff
= static_cast<real
*>(ems
.f
.rvec_array()[0]);
1868 // calculate line gradient in position A
1869 for (gpa
= 0, i
= 0; i
< n
; i
++)
1871 gpa
-= s
[i
] * ff
[i
];
1874 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1875 * relative change in coordinate is smaller than precision
1877 for (minstep
= 0, i
= 0; i
< n
; i
++)
1885 minstep
+= tmp
* tmp
;
1887 minstep
= GMX_REAL_EPS
/ sqrt(minstep
/ n
);
1889 if (stepsize
< minstep
)
1895 // Before taking any steps along the line, store the old position
1897 real
* lastx
= static_cast<real
*>(last
->s
.x
.data()[0]);
1898 real
* lastf
= static_cast<real
*>(last
->f
.data()[0]);
1903 /* Take a step downhill.
1904 * In theory, we should find the actual minimum of the function in this
1905 * direction, somewhere along the line.
1906 * That is quite possible, but it turns out to take 5-10 function evaluations
1907 * for each line. However, we dont really need to find the exact minimum -
1908 * it is much better to start a new BFGS step in a modified direction as soon
1909 * as we are close to it. This will save a lot of energy evaluations.
1911 * In practice, we just try to take a single step.
1912 * If it worked (i.e. lowered the energy), we increase the stepsize but
1913 * continue straight to the next BFGS step without trying to find any minimum,
1914 * i.e. we change the search direction too. If the line was smooth, it is
1915 * likely we are in a smooth region, and then it makes sense to take longer
1916 * steps in the modified search direction too.
1918 * If it didn't work (higher energy), there must be a minimum somewhere between
1919 * the old position and the new one. Then we need to start by finding a lower
1920 * value before we change search direction. Since the energy was apparently
1921 * quite rough, we need to decrease the step size.
1923 * Due to the finite numerical accuracy, it turns out that it is a good idea
1924 * to accept a SMALL increase in energy, if the derivative is still downhill.
1925 * This leads to lower final energies in the tests I've done. / Erik
1928 // State "A" is the first position along the line.
1929 // reference position along line is initially zero
1932 // Check stepsize first. We do not allow displacements
1933 // larger than emstep.
1937 // Pick a new position C by adding stepsize to A.
1940 // Calculate what the largest change in any individual coordinate
1941 // would be (translation along line * gradient along line)
1943 for (i
= 0; i
< n
; i
++)
1946 if (delta
> maxdelta
)
1951 // If any displacement is larger than the stepsize limit, reduce the step
1952 if (maxdelta
> inputrec
->em_stepsize
)
1956 } while (maxdelta
> inputrec
->em_stepsize
);
1958 // Take a trial step and move the coordinate array xc[] to position C
1959 real
* xc
= static_cast<real
*>(sc
->s
.x
.rvec_array()[0]);
1960 for (i
= 0; i
< n
; i
++)
1962 xc
[i
] = lastx
[i
] + c
* s
[i
];
1966 // Calculate energy for the trial step in position C
1967 energyEvaluator
.run(sc
, mu_tot
, vir
, pres
, step
, FALSE
);
1969 // Calc line gradient in position C
1970 real
* fc
= static_cast<real
*>(sc
->f
.rvec_array()[0]);
1971 for (gpc
= 0, i
= 0; i
< n
; i
++)
1973 gpc
-= s
[i
] * fc
[i
]; /* f is negative gradient, thus the sign */
1975 /* Sum the gradient along the line across CPUs */
1978 gmx_sumd(1, &gpc
, cr
);
1981 // This is the max amount of increase in energy we tolerate.
1982 // By allowing VERY small changes (close to numerical precision) we
1983 // frequently find even better (lower) final energies.
1984 tmp
= std::sqrt(GMX_REAL_EPS
) * fabs(sa
->epot
);
1986 // Accept the step if the energy is lower in the new position C (compared to A),
1987 // or if it is not significantly higher and the line derivative is still negative.
1988 foundlower
= sc
->epot
< sa
->epot
|| (gpc
< 0 && sc
->epot
< (sa
->epot
+ tmp
));
1989 // If true, great, we found a better energy. We no longer try to alter the
1990 // stepsize, but simply accept this new better position. The we select a new
1991 // search direction instead, which will be much more efficient than continuing
1992 // to take smaller steps along a line. Set fnorm based on the new C position,
1993 // which will be used to update the stepsize to 1/fnorm further down.
1995 // If false, the energy is NOT lower in point C, i.e. it will be the same
1996 // or higher than in point A. In this case it is pointless to move to point C,
1997 // so we will have to do more iterations along the same line to find a smaller
1998 // value in the interval [A=0.0,C].
1999 // Here, A is still 0.0, but that will change when we do a search in the interval
2000 // [0.0,C] below. That search we will do by interpolation or bisection rather
2001 // than with the stepsize, so no need to modify it. For the next search direction
2002 // it will be reset to 1/fnorm anyway.
2006 // OK, if we didn't find a lower value we will have to locate one now - there must
2007 // be one in the interval [a,c].
2008 // The same thing is valid here, though: Don't spend dozens of iterations to find
2009 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2010 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2011 // I also have a safeguard for potentially really pathological functions so we never
2012 // take more than 20 steps before we give up.
2013 // If we already found a lower value we just skip this step and continue to the update.
2018 // Select a new trial point B in the interval [A,C].
2019 // If the derivatives at points a & c have different sign we interpolate to zero,
2020 // otherwise just do a bisection since there might be multiple minima/maxima
2021 // inside the interval.
2022 if (gpa
< 0 && gpc
> 0)
2024 b
= a
+ gpa
* (a
- c
) / (gpc
- gpa
);
2031 /* safeguard if interpolation close to machine accuracy causes errors:
2032 * never go outside the interval
2034 if (b
<= a
|| b
>= c
)
2039 // Take a trial step to point B
2040 real
* xb
= static_cast<real
*>(sb
->s
.x
.rvec_array()[0]);
2041 for (i
= 0; i
< n
; i
++)
2043 xb
[i
] = lastx
[i
] + b
* s
[i
];
2047 // Calculate energy for the trial step in point B
2048 energyEvaluator
.run(sb
, mu_tot
, vir
, pres
, step
, FALSE
);
2051 // Calculate gradient in point B
2052 real
* fb
= static_cast<real
*>(sb
->f
.rvec_array()[0]);
2053 for (gpb
= 0, i
= 0; i
< n
; i
++)
2055 gpb
-= s
[i
] * fb
[i
]; /* f is negative gradient, thus the sign */
2057 /* Sum the gradient along the line across CPUs */
2060 gmx_sumd(1, &gpb
, cr
);
2063 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2064 // at the new point B, and rename the endpoints of this new interval A and C.
2067 /* Replace c endpoint with b */
2069 /* copy state b to c */
2074 /* Replace a endpoint with b */
2076 /* copy state b to a */
2081 * Stop search as soon as we find a value smaller than the endpoints,
2082 * or if the tolerance is below machine precision.
2083 * Never run more than 20 steps, no matter what.
2086 } while ((sb
->epot
> sa
->epot
|| sb
->epot
> sc
->epot
) && (nminstep
< 20));
2088 if (std::fabs(sb
->epot
- Epot0
) < GMX_REAL_EPS
|| nminstep
>= 20)
2090 /* OK. We couldn't find a significantly lower energy.
2091 * If ncorr==0 this was steepest descent, and then we give up.
2092 * If not, reset memory to restart as steepest descent before quitting.
2104 /* Search in gradient direction */
2105 for (i
= 0; i
< n
; i
++)
2107 dx
[point
][i
] = ff
[i
];
2109 /* Reset stepsize */
2110 stepsize
= 1.0 / fnorm
;
2115 /* Select min energy state of A & C, put the best in xx/ff/Epot
2117 if (sc
->epot
< sa
->epot
)
2138 /* Update the memory information, and calculate a new
2139 * approximation of the inverse hessian
2142 /* Have new data in Epot, xx, ff */
2143 if (ncorr
< nmaxcorr
)
2148 for (i
= 0; i
< n
; i
++)
2150 dg
[point
][i
] = lastf
[i
] - ff
[i
];
2151 dx
[point
][i
] *= step_taken
;
2156 for (i
= 0; i
< n
; i
++)
2158 dgdg
+= dg
[point
][i
] * dg
[point
][i
];
2159 dgdx
+= dg
[point
][i
] * dx
[point
][i
];
2164 rho
[point
] = 1.0 / dgdx
;
2167 if (point
>= nmaxcorr
)
2173 for (i
= 0; i
< n
; i
++)
2180 /* Recursive update. First go back over the memory points */
2181 for (k
= 0; k
< ncorr
; k
++)
2190 for (i
= 0; i
< n
; i
++)
2192 sq
+= dx
[cp
][i
] * p
[i
];
2195 alpha
[cp
] = rho
[cp
] * sq
;
2197 for (i
= 0; i
< n
; i
++)
2199 p
[i
] -= alpha
[cp
] * dg
[cp
][i
];
2203 for (i
= 0; i
< n
; i
++)
2208 /* And then go forward again */
2209 for (k
= 0; k
< ncorr
; k
++)
2212 for (i
= 0; i
< n
; i
++)
2214 yr
+= p
[i
] * dg
[cp
][i
];
2217 beta
= rho
[cp
] * yr
;
2218 beta
= alpha
[cp
] - beta
;
2220 for (i
= 0; i
< n
; i
++)
2222 p
[i
] += beta
* dx
[cp
][i
];
2232 for (i
= 0; i
< n
; i
++)
2236 dx
[point
][i
] = p
[i
];
2244 /* Print it if necessary */
2247 if (mdrunOptions
.verbose
)
2249 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
2250 fprintf(stderr
, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n", step
,
2251 ems
.epot
, ems
.fnorm
/ sqrtNumAtoms
, ems
.fmax
, ems
.a_fmax
+ 1);
2254 /* Store the new (lower) energies */
2255 matrix nullBox
= {};
2256 energyOutput
.addDataAtEnergyStep(false, false, static_cast<double>(step
), mdatoms
->tmass
,
2257 enerd
, nullptr, nullptr, nullptr, nullBox
, nullptr,
2258 nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
2260 do_log
= do_per_step(step
, inputrec
->nstlog
);
2261 do_ene
= do_per_step(step
, inputrec
->nstenergy
);
2263 imdSession
->fillEnergyRecord(step
, TRUE
);
2267 EnergyOutput::printHeader(fplog
, step
, step
);
2269 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), do_ene
, FALSE
, FALSE
,
2270 do_log
? fplog
: nullptr, step
, step
, fcd
, nullptr);
2273 /* Send x and E to IMD client, if bIMD is TRUE. */
2274 if (imdSession
->run(step
, TRUE
, state_global
->box
, state_global
->x
.rvec_array(), 0) && MASTER(cr
))
2276 imdSession
->sendPositionsAndEnergies();
2279 // Reset stepsize in we are doing more iterations
2282 /* Stop when the maximum force lies below tolerance.
2283 * If we have reached machine precision, converged is already set to true.
2285 converged
= converged
|| (ems
.fmax
< inputrec
->em_tol
);
2287 } /* End of the loop */
2291 step
--; /* we never took that last step in this case */
2293 if (ems
.fmax
> inputrec
->em_tol
)
2297 warn_step(fplog
, inputrec
->em_tol
, ems
.fmax
, step
- 1 == number_steps
, FALSE
);
2302 /* If we printed energy and/or logfile last step (which was the last step)
2303 * we don't have to do it again, but otherwise print the final values.
2305 if (!do_log
) /* Write final value to log since we didn't do anythin last step */
2307 EnergyOutput::printHeader(fplog
, step
, step
);
2309 if (!do_ene
|| !do_log
) /* Write final energy file entries */
2311 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), !do_ene
, FALSE
, FALSE
,
2312 !do_log
? fplog
: nullptr, step
, step
, fcd
, nullptr);
2315 /* Print some stuff... */
2318 fprintf(stderr
, "\nwriting lowest energy coordinates.\n");
2322 * For accurate normal mode calculation it is imperative that we
2323 * store the last conformation into the full precision binary trajectory.
2325 * However, we should only do it if we did NOT already write this step
2326 * above (which we did if do_x or do_f was true).
2328 do_x
= !do_per_step(step
, inputrec
->nstxout
);
2329 do_f
= !do_per_step(step
, inputrec
->nstfout
);
2330 write_em_traj(fplog
, cr
, outf
, do_x
, do_f
, ftp2fn(efSTO
, nfile
, fnm
), top_global
, inputrec
,
2331 step
, &ems
, state_global
, observablesHistory
);
2335 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
2336 print_converged(stderr
, LBFGS
, inputrec
->em_tol
, step
, converged
, number_steps
, &ems
, sqrtNumAtoms
);
2337 print_converged(fplog
, LBFGS
, inputrec
->em_tol
, step
, converged
, number_steps
, &ems
, sqrtNumAtoms
);
2339 fprintf(fplog
, "\nPerformed %d energy evaluations in total.\n", neval
);
2342 finish_em(cr
, outf
, walltime_accounting
, wcycle
);
2344 /* To print the actual number of steps we needed somewhere */
2345 walltime_accounting_set_nsteps_done(walltime_accounting
, step
);
2348 void LegacySimulator::do_steep()
2350 const char* SD
= "Steepest Descents";
2351 gmx_localtop_t
top(top_global
->ffparams
);
2352 gmx_global_stat_t gstat
;
2355 gmx_bool bDone
, bAbort
, do_x
, do_f
;
2357 rvec mu_tot
= { 0 };
2360 int steps_accepted
= 0;
2361 auto mdatoms
= mdAtoms
->mdatoms();
2366 "Note that activating steepest-descent energy minimization via the "
2367 "integrator .mdp option and the command gmx mdrun may "
2368 "be available in a different form in a future version of GROMACS, "
2369 "e.g. gmx minimize and an .mdp option.");
2371 /* Create 2 states on the stack and extract pointers that we will swap */
2372 em_state_t s0
{}, s1
{};
2373 em_state_t
* s_min
= &s0
;
2374 em_state_t
* s_try
= &s1
;
2376 /* Init em and store the local state in s_try */
2377 init_em(fplog
, mdlog
, SD
, cr
, inputrec
, imdSession
, pull_work
, state_global
, top_global
, s_try
,
2378 &top
, nrnb
, fr
, mdAtoms
, &gstat
, vsite
, constr
, nullptr);
2379 const bool simulationsShareState
= false;
2380 gmx_mdoutf
* outf
= init_mdoutf(fplog
, nfile
, fnm
, mdrunOptions
, cr
, outputProvider
,
2381 mdModulesNotifier
, inputrec
, top_global
, nullptr, wcycle
,
2382 StartingBehavior::NewSimulation
, simulationsShareState
, ms
);
2383 gmx::EnergyOutput
energyOutput(mdoutf_get_fp_ene(outf
), top_global
, inputrec
, pull_work
, nullptr,
2384 false, StartingBehavior::NewSimulation
, mdModulesNotifier
);
2386 /* Print to log file */
2387 print_em_start(fplog
, cr
, walltime_accounting
, wcycle
, SD
);
2389 /* Set variables for stepsize (in nm). This is the largest
2390 * step that we are going to make in any direction.
2392 ustep
= inputrec
->em_stepsize
;
2395 /* Max number of steps */
2396 nsteps
= inputrec
->nsteps
;
2400 /* Print to the screen */
2401 sp_header(stderr
, SD
, inputrec
->em_tol
, nsteps
);
2405 sp_header(fplog
, SD
, inputrec
->em_tol
, nsteps
);
2407 EnergyEvaluator energyEvaluator
{
2408 fplog
, mdlog
, cr
, ms
, top_global
, &top
, inputrec
, imdSession
, pull_work
,
2409 nrnb
, wcycle
, gstat
, vsite
, constr
, fcd
, mdAtoms
, fr
, runScheduleWork
,
2413 /**** HERE STARTS THE LOOP ****
2414 * count is the counter for the number of steps
2415 * bDone will be TRUE when the minimization has converged
2416 * bAbort will be TRUE when nsteps steps have been performed or when
2417 * the stepsize becomes smaller than is reasonable for machine precision
2422 while (!bDone
&& !bAbort
)
2424 bAbort
= (nsteps
>= 0) && (count
== nsteps
);
2426 /* set new coordinates, except for first step */
2427 bool validStep
= true;
2430 validStep
= do_em_step(cr
, inputrec
, mdatoms
, s_min
, stepsize
, &s_min
->f
, s_try
, constr
, count
);
2435 energyEvaluator
.run(s_try
, mu_tot
, vir
, pres
, count
, count
== 0);
2439 // Signal constraint error during stepping with energy=inf
2440 s_try
->epot
= std::numeric_limits
<real
>::infinity();
2445 EnergyOutput::printHeader(fplog
, count
, count
);
2450 s_min
->epot
= s_try
->epot
;
2453 /* Print it if necessary */
2456 if (mdrunOptions
.verbose
)
2458 fprintf(stderr
, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2459 count
, ustep
, s_try
->epot
, s_try
->fmax
, s_try
->a_fmax
+ 1,
2460 ((count
== 0) || (s_try
->epot
< s_min
->epot
)) ? '\n' : '\r');
2464 if ((count
== 0) || (s_try
->epot
< s_min
->epot
))
2466 /* Store the new (lower) energies */
2467 matrix nullBox
= {};
2468 energyOutput
.addDataAtEnergyStep(false, false, static_cast<double>(count
), mdatoms
->tmass
,
2469 enerd
, nullptr, nullptr, nullptr, nullBox
, nullptr,
2470 nullptr, vir
, pres
, nullptr, mu_tot
, constr
);
2472 imdSession
->fillEnergyRecord(count
, TRUE
);
2474 const bool do_dr
= do_per_step(steps_accepted
, inputrec
->nstdisreout
);
2475 const bool do_or
= do_per_step(steps_accepted
, inputrec
->nstorireout
);
2476 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), TRUE
, do_dr
, do_or
,
2477 fplog
, count
, count
, fcd
, nullptr);
2482 /* Now if the new energy is smaller than the previous...
2483 * or if this is the first step!
2484 * or if we did random steps!
2487 if ((count
== 0) || (s_try
->epot
< s_min
->epot
))
2491 /* Test whether the convergence criterion is met... */
2492 bDone
= (s_try
->fmax
< inputrec
->em_tol
);
2494 /* Copy the arrays for force, positions and energy */
2495 /* The 'Min' array always holds the coords and forces of the minimal
2497 swap_em_state(&s_min
, &s_try
);
2503 /* Write to trn, if necessary */
2504 do_x
= do_per_step(steps_accepted
, inputrec
->nstxout
);
2505 do_f
= do_per_step(steps_accepted
, inputrec
->nstfout
);
2506 write_em_traj(fplog
, cr
, outf
, do_x
, do_f
, nullptr, top_global
, inputrec
, count
, s_min
,
2507 state_global
, observablesHistory
);
2511 /* If energy is not smaller make the step smaller... */
2514 if (DOMAINDECOMP(cr
) && s_min
->s
.ddp_count
!= cr
->dd
->ddp_count
)
2516 /* Reload the old state */
2517 em_dd_partition_system(fplog
, mdlog
, count
, cr
, top_global
, inputrec
, imdSession
,
2518 pull_work
, s_min
, &top
, mdAtoms
, fr
, vsite
, constr
, nrnb
, wcycle
);
2522 // If the force is very small after finishing minimization,
2523 // we risk dividing by zero when calculating the step size.
2524 // So we check first if the minimization has stopped before
2525 // trying to obtain a new step size.
2528 /* Determine new step */
2529 stepsize
= ustep
/ s_min
->fmax
;
2532 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2534 if (count
== nsteps
|| ustep
< 1e-12)
2536 if (count
== nsteps
|| ustep
< 1e-6)
2541 warn_step(fplog
, inputrec
->em_tol
, s_min
->fmax
, count
== nsteps
, constr
!= nullptr);
2546 /* Send IMD energies and positions, if bIMD is TRUE. */
2547 if (imdSession
->run(count
, TRUE
, state_global
->box
,
2548 MASTER(cr
) ? state_global
->x
.rvec_array() : nullptr, 0)
2551 imdSession
->sendPositionsAndEnergies();
2555 } /* End of the loop */
2557 /* Print some data... */
2560 fprintf(stderr
, "\nwriting lowest energy coordinates.\n");
2562 write_em_traj(fplog
, cr
, outf
, TRUE
, inputrec
->nstfout
!= 0, ftp2fn(efSTO
, nfile
, fnm
),
2563 top_global
, inputrec
, count
, s_min
, state_global
, observablesHistory
);
2567 double sqrtNumAtoms
= sqrt(static_cast<double>(state_global
->natoms
));
2569 print_converged(stderr
, SD
, inputrec
->em_tol
, count
, bDone
, nsteps
, s_min
, sqrtNumAtoms
);
2570 print_converged(fplog
, SD
, inputrec
->em_tol
, count
, bDone
, nsteps
, s_min
, sqrtNumAtoms
);
2573 finish_em(cr
, outf
, walltime_accounting
, wcycle
);
2575 /* To print the actual number of steps we needed somewhere */
2576 inputrec
->nsteps
= count
;
2578 walltime_accounting_set_nsteps_done(walltime_accounting
, count
);
2581 void LegacySimulator::do_nm()
2583 const char* NM
= "Normal Mode Analysis";
2585 gmx_localtop_t
top(top_global
->ffparams
);
2586 gmx_global_stat_t gstat
;
2588 rvec mu_tot
= { 0 };
2590 gmx_bool bSparse
; /* use sparse matrix storage format */
2592 gmx_sparsematrix_t
* sparse_matrix
= nullptr;
2593 real
* full_matrix
= nullptr;
2595 /* added with respect to mdrun */
2597 real der_range
= 10.0 * std::sqrt(GMX_REAL_EPS
);
2599 bool bIsMaster
= MASTER(cr
);
2600 auto mdatoms
= mdAtoms
->mdatoms();
2605 "Note that activating normal-mode analysis via the integrator "
2606 ".mdp option and the command gmx mdrun may "
2607 "be available in a different form in a future version of GROMACS, "
2608 "e.g. gmx normal-modes.");
2610 if (constr
!= nullptr)
2614 "Constraints present with Normal Mode Analysis, this combination is not supported");
2617 gmx_shellfc_t
* shellfc
;
2619 em_state_t state_work
{};
2621 /* Init em and store the local state in state_minimum */
2622 init_em(fplog
, mdlog
, NM
, cr
, inputrec
, imdSession
, pull_work
, state_global
, top_global
,
2623 &state_work
, &top
, nrnb
, fr
, mdAtoms
, &gstat
, vsite
, constr
, &shellfc
);
2624 const bool simulationsShareState
= false;
2625 gmx_mdoutf
* outf
= init_mdoutf(fplog
, nfile
, fnm
, mdrunOptions
, cr
, outputProvider
,
2626 mdModulesNotifier
, inputrec
, top_global
, nullptr, wcycle
,
2627 StartingBehavior::NewSimulation
, simulationsShareState
, ms
);
2629 std::vector
<int> atom_index
= get_atom_index(top_global
);
2630 std::vector
<gmx::RVec
> fneg(atom_index
.size(), { 0, 0, 0 });
2631 snew(dfdx
, atom_index
.size());
2637 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2638 " which MIGHT not be accurate enough for normal mode analysis.\n"
2639 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2640 " are fairly modest even if you recompile in double precision.\n\n");
2644 /* Check if we can/should use sparse storage format.
2646 * Sparse format is only useful when the Hessian itself is sparse, which it
2647 * will be when we use a cutoff.
2648 * For small systems (n<1000) it is easier to always use full matrix format, though.
2650 if (EEL_FULL(fr
->ic
->eeltype
) || fr
->rlist
== 0.0)
2652 GMX_LOG(mdlog
.warning
)
2653 .appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2656 else if (atom_index
.size() < 1000)
2658 GMX_LOG(mdlog
.warning
)
2659 .appendTextFormatted("Small system size (N=%zu), using full Hessian format.",
2665 GMX_LOG(mdlog
.warning
).appendText("Using compressed symmetric sparse Hessian format.");
2669 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2670 sz
= DIM
* atom_index
.size();
2672 fprintf(stderr
, "Allocating Hessian memory...\n\n");
2676 sparse_matrix
= gmx_sparsematrix_init(sz
);
2677 sparse_matrix
->compressed_symmetric
= TRUE
;
2681 snew(full_matrix
, sz
* sz
);
2684 /* Write start time and temperature */
2685 print_em_start(fplog
, cr
, walltime_accounting
, wcycle
, NM
);
2687 /* fudge nr of steps to nr of atoms */
2688 inputrec
->nsteps
= atom_index
.size() * 2;
2692 fprintf(stderr
, "starting normal mode calculation '%s'\n%" PRId64
" steps.\n\n",
2693 *(top_global
->name
), inputrec
->nsteps
);
2696 nnodes
= cr
->nnodes
;
2698 /* Make evaluate_energy do a single node force calculation */
2700 EnergyEvaluator energyEvaluator
{
2701 fplog
, mdlog
, cr
, ms
, top_global
, &top
, inputrec
, imdSession
, pull_work
,
2702 nrnb
, wcycle
, gstat
, vsite
, constr
, fcd
, mdAtoms
, fr
, runScheduleWork
,
2705 energyEvaluator
.run(&state_work
, mu_tot
, vir
, pres
, -1, TRUE
);
2706 cr
->nnodes
= nnodes
;
2708 /* if forces are not small, warn user */
2709 get_state_f_norm_max(cr
, &(inputrec
->opts
), mdatoms
, &state_work
);
2711 GMX_LOG(mdlog
.warning
).appendTextFormatted("Maximum force:%12.5e", state_work
.fmax
);
2712 if (state_work
.fmax
> 1.0e-3)
2714 GMX_LOG(mdlog
.warning
)
2716 "The force is probably not small enough to "
2717 "ensure that you are at a minimum.\n"
2718 "Be aware that negative eigenvalues may occur\n"
2719 "when the resulting matrix is diagonalized.");
2722 /***********************************************************
2724 * Loop over all pairs in matrix
2726 * do_force called twice. Once with positive and
2727 * once with negative displacement
2729 ************************************************************/
2731 /* Steps are divided one by one over the nodes */
2733 auto state_work_x
= makeArrayRef(state_work
.s
.x
);
2734 auto state_work_f
= makeArrayRef(state_work
.f
);
2735 for (index aid
= cr
->nodeid
; aid
< ssize(atom_index
); aid
+= nnodes
)
2737 size_t atom
= atom_index
[aid
];
2738 for (size_t d
= 0; d
< DIM
; d
++)
2741 int force_flags
= GMX_FORCE_STATECHANGED
| GMX_FORCE_ALLFORCES
;
2744 x_min
= state_work_x
[atom
][d
];
2746 for (unsigned int dx
= 0; (dx
< 2); dx
++)
2750 state_work_x
[atom
][d
] = x_min
- der_range
;
2754 state_work_x
[atom
][d
] = x_min
+ der_range
;
2757 /* Make evaluate_energy do a single node force calculation */
2761 /* Now is the time to relax the shells */
2762 relax_shell_flexcon(
2763 fplog
, cr
, ms
, mdrunOptions
.verbose
, nullptr, step
, inputrec
, imdSession
,
2764 pull_work
, bNS
, force_flags
, &top
, constr
, enerd
, fcd
, state_work
.s
.natoms
,
2765 state_work
.s
.x
.arrayRefWithPadding(), state_work
.s
.v
.arrayRefWithPadding(),
2766 state_work
.s
.box
, state_work
.s
.lambda
, &state_work
.s
.hist
,
2767 state_work
.f
.arrayRefWithPadding(), vir
, mdatoms
, nrnb
, wcycle
, shellfc
,
2768 fr
, runScheduleWork
, t
, mu_tot
, vsite
, DDBalanceRegionHandler(nullptr));
2774 energyEvaluator
.run(&state_work
, mu_tot
, vir
, pres
, aid
* 2 + dx
, FALSE
);
2777 cr
->nnodes
= nnodes
;
2781 std::copy(state_work_f
.begin(), state_work_f
.begin() + atom_index
.size(),
2786 /* x is restored to original */
2787 state_work_x
[atom
][d
] = x_min
;
2789 for (size_t j
= 0; j
< atom_index
.size(); j
++)
2791 for (size_t k
= 0; (k
< DIM
); k
++)
2793 dfdx
[j
][k
] = -(state_work_f
[atom_index
[j
]][k
] - fneg
[j
][k
]) / (2 * der_range
);
2800 # define mpi_type GMX_MPI_REAL
2801 MPI_Send(dfdx
[0], atom_index
.size() * DIM
, mpi_type
, MASTER(cr
), cr
->nodeid
,
2802 cr
->mpi_comm_mygroup
);
2807 for (index node
= 0; (node
< nnodes
&& aid
+ node
< ssize(atom_index
)); node
++)
2813 MPI_Recv(dfdx
[0], atom_index
.size() * DIM
, mpi_type
, node
, node
,
2814 cr
->mpi_comm_mygroup
, &stat
);
2819 row
= (aid
+ node
) * DIM
+ d
;
2821 for (size_t j
= 0; j
< atom_index
.size(); j
++)
2823 for (size_t k
= 0; k
< DIM
; k
++)
2829 if (col
>= row
&& dfdx
[j
][k
] != 0.0)
2831 gmx_sparsematrix_increment_value(sparse_matrix
, row
, col
, dfdx
[j
][k
]);
2836 full_matrix
[row
* sz
+ col
] = dfdx
[j
][k
];
2843 if (mdrunOptions
.verbose
&& fplog
)
2848 /* write progress */
2849 if (bIsMaster
&& mdrunOptions
.verbose
)
2851 fprintf(stderr
, "\rFinished step %d out of %td",
2852 std::min
<int>(atom
+ nnodes
, atom_index
.size()), ssize(atom_index
));
2859 fprintf(stderr
, "\n\nWriting Hessian...\n");
2860 gmx_mtxio_write(ftp2fn(efMTX
, nfile
, fnm
), sz
, sz
, full_matrix
, sparse_matrix
);
2863 finish_em(cr
, outf
, walltime_accounting
, wcycle
);
2865 walltime_accounting_set_nsteps_done(walltime_accounting
, atom_index
.size() * 2);