Call Constraints::setConstraints with other setup code
[gromacs.git] / src / gromacs / mdrun / minimize.cpp
blob96f0a6c3c40acae8f413913558ea6a25f1e5b047
1 /*
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, 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.
37 /*! \internal \file
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
45 #include "gmxpre.h"
47 #include "config.h"
49 #include <cmath>
50 #include <cstring>
51 #include <ctime>
53 #include <algorithm>
54 #include <vector>
56 #include "gromacs/commandline/filenm.h"
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/ewald/pme.h"
60 #include "gromacs/fileio/confio.h"
61 #include "gromacs/fileio/mtxio.h"
62 #include "gromacs/gmxlib/network.h"
63 #include "gromacs/gmxlib/nrnb.h"
64 #include "gromacs/imd/imd.h"
65 #include "gromacs/linearalgebra/sparsematrix.h"
66 #include "gromacs/listed-forces/manage-threading.h"
67 #include "gromacs/math/functions.h"
68 #include "gromacs/math/vec.h"
69 #include "gromacs/mdlib/constr.h"
70 #include "gromacs/mdlib/force.h"
71 #include "gromacs/mdlib/forcerec.h"
72 #include "gromacs/mdlib/gmx_omp_nthreads.h"
73 #include "gromacs/mdlib/md_support.h"
74 #include "gromacs/mdlib/mdatoms.h"
75 #include "gromacs/mdlib/mdebin.h"
76 #include "gromacs/mdlib/mdrun.h"
77 #include "gromacs/mdlib/mdsetup.h"
78 #include "gromacs/mdlib/ns.h"
79 #include "gromacs/mdlib/shellfc.h"
80 #include "gromacs/mdlib/sim_util.h"
81 #include "gromacs/mdlib/tgroup.h"
82 #include "gromacs/mdlib/trajectory_writing.h"
83 #include "gromacs/mdlib/update.h"
84 #include "gromacs/mdlib/vsite.h"
85 #include "gromacs/mdtypes/commrec.h"
86 #include "gromacs/mdtypes/inputrec.h"
87 #include "gromacs/mdtypes/md_enums.h"
88 #include "gromacs/mdtypes/state.h"
89 #include "gromacs/pbcutil/mshift.h"
90 #include "gromacs/pbcutil/pbc.h"
91 #include "gromacs/timing/wallcycle.h"
92 #include "gromacs/timing/walltime_accounting.h"
93 #include "gromacs/topology/mtop_util.h"
94 #include "gromacs/topology/topology.h"
95 #include "gromacs/utility/cstringutil.h"
96 #include "gromacs/utility/exceptions.h"
97 #include "gromacs/utility/fatalerror.h"
98 #include "gromacs/utility/logger.h"
99 #include "gromacs/utility/smalloc.h"
101 #include "integrator.h"
103 //! Utility structure for manipulating states during EM
104 typedef struct {
105 //! Copy of the global state
106 t_state s;
107 //! Force array
108 PaddedRVecVector f;
109 //! Potential energy
110 real epot;
111 //! Norm of the force
112 real fnorm;
113 //! Maximum force
114 real fmax;
115 //! Direction
116 int a_fmax;
117 } em_state_t;
119 //! Print the EM starting conditions
120 static void print_em_start(FILE *fplog,
121 const t_commrec *cr,
122 gmx_walltime_accounting_t walltime_accounting,
123 gmx_wallcycle_t wcycle,
124 const char *name)
126 walltime_accounting_start(walltime_accounting);
127 wallcycle_start(wcycle, ewcRUN);
128 print_start(fplog, cr, walltime_accounting, name);
131 //! Stop counting time for EM
132 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
133 gmx_wallcycle_t wcycle)
135 wallcycle_stop(wcycle, ewcRUN);
137 walltime_accounting_end(walltime_accounting);
140 //! Printing a log file and console header
141 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
143 fprintf(out, "\n");
144 fprintf(out, "%s:\n", minimizer);
145 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
146 fprintf(out, " Number of steps = %12d\n", nsteps);
149 //! Print warning message
150 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
152 char buffer[2048];
153 if (bLastStep)
155 sprintf(buffer,
156 "\nEnergy minimization reached the maximum number "
157 "of steps before the forces reached the requested "
158 "precision Fmax < %g.\n", ftol);
160 else
162 sprintf(buffer,
163 "\nEnergy minimization has stopped, but the forces have "
164 "not converged to the requested precision Fmax < %g (which "
165 "may not be possible for your system). It stopped "
166 "because the algorithm tried to make a new step whose size "
167 "was too small, or there was no change in the energy since "
168 "last step. Either way, we regard the minimization as "
169 "converged to within the available machine precision, "
170 "given your starting configuration and EM parameters.\n%s%s",
171 ftol,
172 sizeof(real) < sizeof(double) ?
173 "\nDouble precision normally gives you higher accuracy, but "
174 "this is often not needed for preparing to run molecular "
175 "dynamics.\n" :
177 bConstrain ?
178 "You might need to increase your constraint accuracy, or turn\n"
179 "off constraints altogether (set constraints = none in mdp file)\n" :
180 "");
182 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
185 //! Print message about convergence of the EM
186 static void print_converged(FILE *fp, const char *alg, real ftol,
187 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
188 const em_state_t *ems, double sqrtNumAtoms)
190 char buf[STEPSTRSIZE];
192 if (bDone)
194 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
195 alg, ftol, gmx_step_str(count, buf));
197 else if (count < nsteps)
199 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
200 "but did not reach the requested Fmax < %g.\n",
201 alg, gmx_step_str(count, buf), ftol);
203 else
205 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
206 alg, ftol, gmx_step_str(count, buf));
209 #if GMX_DOUBLE
210 fprintf(fp, "Potential Energy = %21.14e\n", ems->epot);
211 fprintf(fp, "Maximum force = %21.14e on atom %d\n", ems->fmax, ems->a_fmax + 1);
212 fprintf(fp, "Norm of force = %21.14e\n", ems->fnorm/sqrtNumAtoms);
213 #else
214 fprintf(fp, "Potential Energy = %14.7e\n", ems->epot);
215 fprintf(fp, "Maximum force = %14.7e on atom %d\n", ems->fmax, ems->a_fmax + 1);
216 fprintf(fp, "Norm of force = %14.7e\n", ems->fnorm/sqrtNumAtoms);
217 #endif
220 //! Compute the norm and max of the force array in parallel
221 static void get_f_norm_max(const t_commrec *cr,
222 t_grpopts *opts, t_mdatoms *mdatoms, const rvec *f,
223 real *fnorm, real *fmax, int *a_fmax)
225 double fnorm2, *sum;
226 real fmax2, fam;
227 int la_max, a_max, start, end, i, m, gf;
229 /* This routine finds the largest force and returns it.
230 * On parallel machines the global max is taken.
232 fnorm2 = 0;
233 fmax2 = 0;
234 la_max = -1;
235 start = 0;
236 end = mdatoms->homenr;
237 if (mdatoms->cFREEZE)
239 for (i = start; i < end; i++)
241 gf = mdatoms->cFREEZE[i];
242 fam = 0;
243 for (m = 0; m < DIM; m++)
245 if (!opts->nFreeze[gf][m])
247 fam += gmx::square(f[i][m]);
250 fnorm2 += fam;
251 if (fam > fmax2)
253 fmax2 = fam;
254 la_max = i;
258 else
260 for (i = start; i < end; i++)
262 fam = norm2(f[i]);
263 fnorm2 += fam;
264 if (fam > fmax2)
266 fmax2 = fam;
267 la_max = i;
272 if (la_max >= 0 && DOMAINDECOMP(cr))
274 a_max = cr->dd->globalAtomIndices[la_max];
276 else
278 a_max = la_max;
280 if (PAR(cr))
282 snew(sum, 2*cr->nnodes+1);
283 sum[2*cr->nodeid] = fmax2;
284 sum[2*cr->nodeid+1] = a_max;
285 sum[2*cr->nnodes] = fnorm2;
286 gmx_sumd(2*cr->nnodes+1, sum, cr);
287 fnorm2 = sum[2*cr->nnodes];
288 /* Determine the global maximum */
289 for (i = 0; i < cr->nnodes; i++)
291 if (sum[2*i] > fmax2)
293 fmax2 = sum[2*i];
294 a_max = (int)(sum[2*i+1] + 0.5);
297 sfree(sum);
300 if (fnorm)
302 *fnorm = sqrt(fnorm2);
304 if (fmax)
306 *fmax = sqrt(fmax2);
308 if (a_fmax)
310 *a_fmax = a_max;
314 //! Compute the norm of the force
315 static void get_state_f_norm_max(const t_commrec *cr,
316 t_grpopts *opts, t_mdatoms *mdatoms,
317 em_state_t *ems)
319 get_f_norm_max(cr, opts, mdatoms, as_rvec_array(ems->f.data()),
320 &ems->fnorm, &ems->fmax, &ems->a_fmax);
323 //! Initialize the energy minimization
324 static void init_em(FILE *fplog, const char *title,
325 const t_commrec *cr,
326 const gmx_multisim_t *ms,
327 gmx::IMDOutputProvider *outputProvider,
328 t_inputrec *ir,
329 const MdrunOptions &mdrunOptions,
330 t_state *state_global, gmx_mtop_t *top_global,
331 em_state_t *ems, gmx_localtop_t **top,
332 t_nrnb *nrnb, rvec mu_tot,
333 t_forcerec *fr, gmx_enerdata_t **enerd,
334 t_graph **graph, gmx::MDAtoms *mdAtoms, gmx_global_stat_t *gstat,
335 gmx_vsite_t *vsite, gmx::Constraints *constr, gmx_shellfc_t **shellfc,
336 int nfile, const t_filenm fnm[],
337 gmx_mdoutf_t *outf, t_mdebin **mdebin,
338 gmx_wallcycle_t wcycle)
340 real dvdl_constr;
342 if (fplog)
344 fprintf(fplog, "Initiating %s\n", title);
347 if (MASTER(cr))
349 state_global->ngtc = 0;
351 /* Initialize lambda variables */
352 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, nullptr);
355 init_nrnb(nrnb);
357 /* Interactive molecular dynamics */
358 init_IMD(ir, cr, ms, top_global, fplog, 1,
359 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
360 nfile, fnm, nullptr, mdrunOptions);
362 if (ir->eI == eiNM)
364 GMX_ASSERT(shellfc != nullptr, "With NM we always support shells");
366 *shellfc = init_shell_flexcon(stdout,
367 top_global,
368 constr ? constr->numFlexibleConstraints() : 0,
369 ir->nstcalcenergy,
370 DOMAINDECOMP(cr));
372 else
374 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
376 /* With energy minimization, shells and flexible constraints are
377 * automatically minimized when treated like normal DOFS.
379 if (shellfc != nullptr)
381 *shellfc = nullptr;
385 auto mdatoms = mdAtoms->mdatoms();
386 if (DOMAINDECOMP(cr))
388 *top = dd_init_local_top(top_global);
390 dd_init_local_state(cr->dd, state_global, &ems->s);
392 /* Distribute the charge groups over the nodes from the master node */
393 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
394 state_global, top_global, ir,
395 &ems->s, &ems->f, mdAtoms, *top,
396 fr, vsite, constr,
397 nrnb, nullptr, FALSE);
398 dd_store_state(cr->dd, &ems->s);
400 *graph = nullptr;
402 else
404 state_change_natoms(state_global, state_global->natoms);
405 /* Just copy the state */
406 ems->s = *state_global;
407 state_change_natoms(&ems->s, ems->s.natoms);
408 /* We need to allocate one element extra, since we might use
409 * (unaligned) 4-wide SIMD loads to access rvec entries.
411 ems->f.resize(gmx::paddedRVecVectorSize(ems->s.natoms));
413 snew(*top, 1);
414 mdAlgorithmsSetupAtomData(cr, ir, top_global, *top, fr,
415 graph, mdAtoms,
416 constr, vsite, shellfc ? *shellfc : nullptr);
418 if (vsite)
420 set_vsite_top(vsite, *top, mdatoms);
424 update_mdatoms(mdAtoms->mdatoms(), ems->s.lambda[efptMASS]);
426 if (constr)
428 // TODO how should this cross-module support dependency be managed?
429 if (ir->eConstrAlg == econtSHAKE &&
430 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
432 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
433 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
436 if (!ir->bContinuation)
438 /* Constrain the starting coordinates */
439 dvdl_constr = 0;
440 constr->apply(TRUE, TRUE,
441 -1, 0, 1.0,
442 as_rvec_array(ems->s.x.data()),
443 as_rvec_array(ems->s.x.data()),
444 nullptr,
445 ems->s.box,
446 ems->s.lambda[efptFEP], &dvdl_constr,
447 nullptr, nullptr, gmx::ConstraintVariable::Positions);
451 if (PAR(cr))
453 *gstat = global_stat_init(ir);
455 else
457 *gstat = nullptr;
460 *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, nullptr, wcycle);
462 snew(*enerd, 1);
463 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
464 *enerd);
466 if (mdebin != nullptr)
468 /* Init bin for energy stuff */
469 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, nullptr);
472 clear_rvec(mu_tot);
473 calc_shifts(ems->s.box, fr->shift_vec);
476 //! Finalize the minimization
477 static void finish_em(const t_commrec *cr, gmx_mdoutf_t outf,
478 gmx_walltime_accounting_t walltime_accounting,
479 gmx_wallcycle_t wcycle)
481 if (!thisRankHasDuty(cr, DUTY_PME))
483 /* Tell the PME only node to finish */
484 gmx_pme_send_finish(cr);
487 done_mdoutf(outf);
489 em_time_end(walltime_accounting, wcycle);
492 //! Swap two different EM states during minimization
493 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
495 em_state_t *tmp;
497 tmp = *ems1;
498 *ems1 = *ems2;
499 *ems2 = tmp;
502 //! Save the EM trajectory
503 static void write_em_traj(FILE *fplog, const t_commrec *cr,
504 gmx_mdoutf_t outf,
505 gmx_bool bX, gmx_bool bF, const char *confout,
506 gmx_mtop_t *top_global,
507 t_inputrec *ir, gmx_int64_t step,
508 em_state_t *state,
509 t_state *state_global,
510 ObservablesHistory *observablesHistory)
512 int mdof_flags = 0;
514 if (bX)
516 mdof_flags |= MDOF_X;
518 if (bF)
520 mdof_flags |= MDOF_F;
523 /* If we want IMD output, set appropriate MDOF flag */
524 if (ir->bIMD)
526 mdof_flags |= MDOF_IMD;
529 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
530 top_global, step, (double)step,
531 &state->s, state_global, observablesHistory,
532 state->f);
534 if (confout != nullptr && MASTER(cr))
536 GMX_RELEASE_ASSERT(bX, "The code below assumes that (with domain decomposition), x is collected to state_global in the call above.");
537 /* With domain decomposition the call above collected the state->s.x
538 * into state_global->x. Without DD we copy the local state pointer.
540 if (!DOMAINDECOMP(cr))
542 state_global = &state->s;
545 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
547 /* Make molecules whole only for confout writing */
548 do_pbc_mtop(fplog, ir->ePBC, state->s.box, top_global,
549 as_rvec_array(state_global->x.data()));
552 write_sto_conf_mtop(confout,
553 *top_global->name, top_global,
554 as_rvec_array(state_global->x.data()), nullptr, ir->ePBC, state->s.box);
558 //! \brief Do one minimization step
560 // \returns true when the step succeeded, false when a constraint error occurred
561 static bool do_em_step(const t_commrec *cr,
562 t_inputrec *ir, t_mdatoms *md,
563 em_state_t *ems1, real a, const PaddedRVecVector *force,
564 em_state_t *ems2,
565 gmx::Constraints *constr,
566 gmx_int64_t count)
569 t_state *s1, *s2;
570 int start, end;
571 real dvdl_constr;
572 int nthreads gmx_unused;
574 bool validStep = true;
576 s1 = &ems1->s;
577 s2 = &ems2->s;
579 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
581 gmx_incons("state mismatch in do_em_step");
584 s2->flags = s1->flags;
586 if (s2->natoms != s1->natoms)
588 state_change_natoms(s2, s1->natoms);
589 /* We need to allocate one element extra, since we might use
590 * (unaligned) 4-wide SIMD loads to access rvec entries.
592 ems2->f.resize(gmx::paddedRVecVectorSize(s2->natoms));
594 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
596 s2->cg_gl.resize(s1->cg_gl.size());
599 copy_mat(s1->box, s2->box);
600 /* Copy free energy state */
601 s2->lambda = s1->lambda;
602 copy_mat(s1->box, s2->box);
604 start = 0;
605 end = md->homenr;
607 // cppcheck-suppress unreadVariable
608 nthreads = gmx_omp_nthreads_get(emntUpdate);
609 #pragma omp parallel num_threads(nthreads)
611 const rvec *x1 = as_rvec_array(s1->x.data());
612 rvec *x2 = as_rvec_array(s2->x.data());
613 const rvec *f = as_rvec_array(force->data());
615 int gf = 0;
616 #pragma omp for schedule(static) nowait
617 for (int i = start; i < end; i++)
621 if (md->cFREEZE)
623 gf = md->cFREEZE[i];
625 for (int m = 0; m < DIM; m++)
627 if (ir->opts.nFreeze[gf][m])
629 x2[i][m] = x1[i][m];
631 else
633 x2[i][m] = x1[i][m] + a*f[i][m];
637 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
640 if (s2->flags & (1<<estCGP))
642 /* Copy the CG p vector */
643 const rvec *p1 = as_rvec_array(s1->cg_p.data());
644 rvec *p2 = as_rvec_array(s2->cg_p.data());
645 #pragma omp for schedule(static) nowait
646 for (int i = start; i < end; i++)
648 // Trivial OpenMP block that does not throw
649 copy_rvec(p1[i], p2[i]);
653 if (DOMAINDECOMP(cr))
655 s2->ddp_count = s1->ddp_count;
657 /* OpenMP does not supported unsigned loop variables */
658 #pragma omp for schedule(static) nowait
659 for (int i = 0; i < static_cast<int>(s2->cg_gl.size()); i++)
661 s2->cg_gl[i] = s1->cg_gl[i];
663 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
667 if (constr)
669 dvdl_constr = 0;
670 validStep =
671 constr->apply(TRUE, TRUE,
672 count, 0, 1.0,
673 as_rvec_array(s1->x.data()), as_rvec_array(s2->x.data()),
674 nullptr, s2->box,
675 s2->lambda[efptBONDED], &dvdl_constr,
676 nullptr, nullptr, gmx::ConstraintVariable::Positions);
678 // We should move this check to the different minimizers
679 if (!validStep && ir->eI != eiSteep)
681 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
682 EI(ir->eI), EI(eiSteep), EI(ir->eI));
686 return validStep;
689 //! Prepare EM for using domain decomposition parallellization
690 static void em_dd_partition_system(FILE *fplog, int step, const t_commrec *cr,
691 gmx_mtop_t *top_global, t_inputrec *ir,
692 em_state_t *ems, gmx_localtop_t *top,
693 gmx::MDAtoms *mdAtoms, t_forcerec *fr,
694 gmx_vsite_t *vsite, gmx::Constraints *constr,
695 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
697 /* Repartition the domain decomposition */
698 dd_partition_system(fplog, step, cr, FALSE, 1,
699 nullptr, top_global, ir,
700 &ems->s, &ems->f,
701 mdAtoms, top, fr, vsite, constr,
702 nrnb, wcycle, FALSE);
703 dd_store_state(cr->dd, &ems->s);
706 namespace
709 /*! \brief Class to handle the work of setting and doing an energy evaluation.
711 * This class is a mere aggregate of parameters to pass to evaluate an
712 * energy, so that future changes to names and types of them consume
713 * less time when refactoring other code.
715 * Aggregate initialization is used, for which the chief risk is that
716 * if a member is added at the end and not all initializer lists are
717 * updated, then the member will be value initialized, which will
718 * typically mean initialization to zero.
720 * We only want to construct one of these with an initializer list, so
721 * we explicitly delete the default constructor. */
722 class EnergyEvaluator
724 public:
725 //! We only intend to construct such objects with an initializer list.
726 #if __GNUC__ > 4 || (__GNUC__ == 4 && __GNUC_MINOR__ >= 9)
727 // Aspects of the C++11 spec changed after GCC 4.8.5, and
728 // compilation of the initializer list construction in
729 // runner.cpp fails in GCC 4.8.5.
730 EnergyEvaluator() = delete;
731 #endif
732 /*! \brief Evaluates an energy on the state in \c ems.
734 * \todo In practice, the same objects mu_tot, vir, and pres
735 * are always passed to this function, so we would rather have
736 * them as data members. However, their C-array types are
737 * unsuited for aggregate initialization. When the types
738 * improve, the call signature of this method can be reduced.
740 void run(em_state_t *ems, rvec mu_tot,
741 tensor vir, tensor pres,
742 gmx_int64_t count, gmx_bool bFirst);
743 //! Handles logging.
744 FILE *fplog;
745 //! Handles communication.
746 const t_commrec *cr;
747 //! Coordinates multi-simulations.
748 const gmx_multisim_t *ms;
749 //! Holds the simulation topology.
750 gmx_mtop_t *top_global;
751 //! Holds the domain topology.
752 gmx_localtop_t *top;
753 //! User input options.
754 t_inputrec *inputrec;
755 //! Manages flop accounting.
756 t_nrnb *nrnb;
757 //! Manages wall cycle accounting.
758 gmx_wallcycle_t wcycle;
759 //! Coordinates global reduction.
760 gmx_global_stat_t gstat;
761 //! Handles virtual sites.
762 gmx_vsite_t *vsite;
763 //! Handles constraints.
764 gmx::Constraints *constr;
765 //! Handles strange things.
766 t_fcdata *fcd;
767 //! Molecular graph for SHAKE.
768 t_graph *graph;
769 //! Per-atom data for this domain.
770 gmx::MDAtoms *mdAtoms;
771 //! Handles how to calculate the forces.
772 t_forcerec *fr;
773 //! Stores the computed energies.
774 gmx_enerdata_t *enerd;
777 void
778 EnergyEvaluator::run(em_state_t *ems, rvec mu_tot,
779 tensor vir, tensor pres,
780 gmx_int64_t count, gmx_bool bFirst)
782 real t;
783 gmx_bool bNS;
784 tensor force_vir, shake_vir, ekin;
785 real dvdl_constr, prescorr, enercorr, dvdlcorr;
786 real terminate = 0;
788 /* Set the time to the initial time, the time does not change during EM */
789 t = inputrec->init_t;
791 if (bFirst ||
792 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
794 /* This is the first state or an old state used before the last ns */
795 bNS = TRUE;
797 else
799 bNS = FALSE;
800 if (inputrec->nstlist > 0)
802 bNS = TRUE;
806 if (vsite)
808 construct_vsites(vsite, as_rvec_array(ems->s.x.data()), 1, nullptr,
809 top->idef.iparams, top->idef.il,
810 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
813 if (DOMAINDECOMP(cr) && bNS)
815 /* Repartition the domain decomposition */
816 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
817 ems, top, mdAtoms, fr, vsite, constr,
818 nrnb, wcycle);
821 /* Calc force & energy on new trial position */
822 /* do_force always puts the charge groups in the box and shifts again
823 * We do not unshift, so molecules are always whole in congrad.c
825 do_force(fplog, cr, ms, inputrec, nullptr,
826 count, nrnb, wcycle, top, &top_global->groups,
827 ems->s.box, ems->s.x, &ems->s.hist,
828 ems->f, force_vir, mdAtoms->mdatoms(), enerd, fcd,
829 ems->s.lambda, graph, fr, vsite, mu_tot, t, nullptr,
830 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
831 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
832 (bNS ? GMX_FORCE_NS : 0),
833 DOMAINDECOMP(cr) ?
834 DdOpenBalanceRegionBeforeForceComputation::yes :
835 DdOpenBalanceRegionBeforeForceComputation::no,
836 DOMAINDECOMP(cr) ?
837 DdCloseBalanceRegionAfterForceComputation::yes :
838 DdCloseBalanceRegionAfterForceComputation::no);
840 /* Clear the unused shake virial and pressure */
841 clear_mat(shake_vir);
842 clear_mat(pres);
844 /* Communicate stuff when parallel */
845 if (PAR(cr) && inputrec->eI != eiNM)
847 wallcycle_start(wcycle, ewcMoveE);
849 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
850 inputrec, nullptr, nullptr, nullptr, 1, &terminate,
851 nullptr, FALSE,
852 CGLO_ENERGY |
853 CGLO_PRESSURE |
854 CGLO_CONSTRAINT);
856 wallcycle_stop(wcycle, ewcMoveE);
859 /* Calculate long range corrections to pressure and energy */
860 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
861 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
862 enerd->term[F_DISPCORR] = enercorr;
863 enerd->term[F_EPOT] += enercorr;
864 enerd->term[F_PRES] += prescorr;
865 enerd->term[F_DVDL] += dvdlcorr;
867 ems->epot = enerd->term[F_EPOT];
869 if (constr)
871 /* Project out the constraint components of the force */
872 dvdl_constr = 0;
873 rvec *f_rvec = as_rvec_array(ems->f.data());
874 constr->apply(FALSE, FALSE,
875 count, 0, 1.0,
876 as_rvec_array(ems->s.x.data()), f_rvec, f_rvec,
877 ems->s.box,
878 ems->s.lambda[efptBONDED], &dvdl_constr,
879 nullptr, &shake_vir, gmx::ConstraintVariable::ForceDispl);
880 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
881 m_add(force_vir, shake_vir, vir);
883 else
885 copy_mat(force_vir, vir);
888 clear_mat(ekin);
889 enerd->term[F_PRES] =
890 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
892 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
894 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
896 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
900 } // namespace
902 //! Parallel utility summing energies and forces
903 static double reorder_partsum(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
904 gmx_mtop_t *top_global,
905 em_state_t *s_min, em_state_t *s_b)
907 t_block *cgs_gl;
908 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
909 double partsum;
910 unsigned char *grpnrFREEZE;
912 if (debug)
914 fprintf(debug, "Doing reorder_partsum\n");
917 const rvec *fm = as_rvec_array(s_min->f.data());
918 const rvec *fb = as_rvec_array(s_b->f.data());
920 cgs_gl = dd_charge_groups_global(cr->dd);
921 index = cgs_gl->index;
923 /* Collect fm in a global vector fmg.
924 * This conflicts with the spirit of domain decomposition,
925 * but to fully optimize this a much more complicated algorithm is required.
927 rvec *fmg;
928 snew(fmg, top_global->natoms);
930 ncg = s_min->s.cg_gl.size();
931 cg_gl = s_min->s.cg_gl.data();
932 i = 0;
933 for (c = 0; c < ncg; c++)
935 cg = cg_gl[c];
936 a0 = index[cg];
937 a1 = index[cg+1];
938 for (a = a0; a < a1; a++)
940 copy_rvec(fm[i], fmg[a]);
941 i++;
944 gmx_sum(top_global->natoms*3, fmg[0], cr);
946 /* Now we will determine the part of the sum for the cgs in state s_b */
947 ncg = s_b->s.cg_gl.size();
948 cg_gl = s_b->s.cg_gl.data();
949 partsum = 0;
950 i = 0;
951 gf = 0;
952 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
953 for (c = 0; c < ncg; c++)
955 cg = cg_gl[c];
956 a0 = index[cg];
957 a1 = index[cg+1];
958 for (a = a0; a < a1; a++)
960 if (mdatoms->cFREEZE && grpnrFREEZE)
962 gf = grpnrFREEZE[i];
964 for (m = 0; m < DIM; m++)
966 if (!opts->nFreeze[gf][m])
968 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
971 i++;
975 sfree(fmg);
977 return partsum;
980 //! Print some stuff, like beta, whatever that means.
981 static real pr_beta(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
982 gmx_mtop_t *top_global,
983 em_state_t *s_min, em_state_t *s_b)
985 double sum;
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 &&
994 s_b->s.ddp_count == cr->dd->ddp_count))
996 const rvec *fm = as_rvec_array(s_min->f.data());
997 const rvec *fb = as_rvec_array(s_b->f.data());
998 sum = 0;
999 int gf = 0;
1000 /* This part of code can be incorrect with DD,
1001 * since the atom ordering in s_b and s_min might differ.
1003 for (int i = 0; i < mdatoms->homenr; i++)
1005 if (mdatoms->cFREEZE)
1007 gf = mdatoms->cFREEZE[i];
1009 for (int m = 0; m < DIM; m++)
1011 if (!opts->nFreeze[gf][m])
1013 sum += (fb[i][m] - fm[i][m])*fb[i][m];
1018 else
1020 /* We need to reorder cgs while summing */
1021 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
1023 if (PAR(cr))
1025 gmx_sumd(1, &sum, cr);
1028 return sum/gmx::square(s_min->fnorm);
1031 namespace gmx
1034 void
1035 Integrator::do_cg()
1037 const char *CG = "Polak-Ribiere Conjugate Gradients";
1039 gmx_localtop_t *top;
1040 gmx_enerdata_t *enerd;
1041 gmx_global_stat_t gstat;
1042 t_graph *graph;
1043 double tmp, minstep;
1044 real stepsize;
1045 real a, b, c, beta = 0.0;
1046 real epot_repl = 0;
1047 real pnorm;
1048 t_mdebin *mdebin;
1049 gmx_bool converged, foundlower;
1050 rvec mu_tot;
1051 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1052 tensor vir, pres;
1053 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1054 gmx_mdoutf_t outf;
1055 int m, step, nminstep;
1056 auto mdatoms = mdAtoms->mdatoms();
1058 step = 0;
1060 // Ensure the extra per-atom state array gets allocated
1061 state_global->flags |= (1<<estCGP);
1063 /* Create 4 states on the stack and extract pointers that we will swap */
1064 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1065 em_state_t *s_min = &s0;
1066 em_state_t *s_a = &s1;
1067 em_state_t *s_b = &s2;
1068 em_state_t *s_c = &s3;
1070 /* Init em and store the local state in s_min */
1071 init_em(fplog, CG, cr, ms, outputProvider, inputrec, mdrunOptions,
1072 state_global, top_global, s_min, &top,
1073 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1074 vsite, constr, nullptr,
1075 nfile, fnm, &outf, &mdebin, wcycle);
1077 /* Print to log file */
1078 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1080 /* Max number of steps */
1081 number_steps = inputrec->nsteps;
1083 if (MASTER(cr))
1085 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1087 if (fplog)
1089 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1092 EnergyEvaluator energyEvaluator {
1093 fplog, cr, ms,
1094 top_global, top,
1095 inputrec, nrnb, wcycle, gstat,
1096 vsite, constr, fcd, graph,
1097 mdAtoms, fr, enerd
1099 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1100 /* do_force always puts the charge groups in the box and shifts again
1101 * We do not unshift, so molecules are always whole in congrad.c
1103 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1105 if (MASTER(cr))
1107 /* Copy stuff to the energy bin for easy printing etc. */
1108 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1109 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1110 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1112 print_ebin_header(fplog, step, step);
1113 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1114 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1117 /* Estimate/guess the initial stepsize */
1118 stepsize = inputrec->em_stepsize/s_min->fnorm;
1120 if (MASTER(cr))
1122 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1123 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1124 s_min->fmax, s_min->a_fmax+1);
1125 fprintf(stderr, " F-Norm = %12.5e\n",
1126 s_min->fnorm/sqrtNumAtoms);
1127 fprintf(stderr, "\n");
1128 /* and copy to the log file too... */
1129 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1130 s_min->fmax, s_min->a_fmax+1);
1131 fprintf(fplog, " F-Norm = %12.5e\n",
1132 s_min->fnorm/sqrtNumAtoms);
1133 fprintf(fplog, "\n");
1135 /* Start the loop over CG steps.
1136 * Each successful step is counted, and we continue until
1137 * we either converge or reach the max number of steps.
1139 converged = FALSE;
1140 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1143 /* start taking steps in a new direction
1144 * First time we enter the routine, beta=0, and the direction is
1145 * simply the negative gradient.
1148 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1149 rvec *pm = as_rvec_array(s_min->s.cg_p.data());
1150 const rvec *sfm = as_rvec_array(s_min->f.data());
1151 double gpa = 0;
1152 int gf = 0;
1153 for (int i = 0; i < mdatoms->homenr; i++)
1155 if (mdatoms->cFREEZE)
1157 gf = mdatoms->cFREEZE[i];
1159 for (m = 0; m < DIM; m++)
1161 if (!inputrec->opts.nFreeze[gf][m])
1163 pm[i][m] = sfm[i][m] + beta*pm[i][m];
1164 gpa -= pm[i][m]*sfm[i][m];
1165 /* f is negative gradient, thus the sign */
1167 else
1169 pm[i][m] = 0;
1174 /* Sum the gradient along the line across CPUs */
1175 if (PAR(cr))
1177 gmx_sumd(1, &gpa, cr);
1180 /* Calculate the norm of the search vector */
1181 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1183 /* Just in case stepsize reaches zero due to numerical precision... */
1184 if (stepsize <= 0)
1186 stepsize = inputrec->em_stepsize/pnorm;
1190 * Double check the value of the derivative in the search direction.
1191 * If it is positive it must be due to the old information in the
1192 * CG formula, so just remove that and start over with beta=0.
1193 * This corresponds to a steepest descent step.
1195 if (gpa > 0)
1197 beta = 0;
1198 step--; /* Don't count this step since we are restarting */
1199 continue; /* Go back to the beginning of the big for-loop */
1202 /* Calculate minimum allowed stepsize, before the average (norm)
1203 * relative change in coordinate is smaller than precision
1205 minstep = 0;
1206 for (int i = 0; i < mdatoms->homenr; i++)
1208 for (m = 0; m < DIM; m++)
1210 tmp = fabs(s_min->s.x[i][m]);
1211 if (tmp < 1.0)
1213 tmp = 1.0;
1215 tmp = pm[i][m]/tmp;
1216 minstep += tmp*tmp;
1219 /* Add up from all CPUs */
1220 if (PAR(cr))
1222 gmx_sumd(1, &minstep, cr);
1225 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1227 if (stepsize < minstep)
1229 converged = TRUE;
1230 break;
1233 /* Write coordinates if necessary */
1234 do_x = do_per_step(step, inputrec->nstxout);
1235 do_f = do_per_step(step, inputrec->nstfout);
1237 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1238 top_global, inputrec, step,
1239 s_min, state_global, observablesHistory);
1241 /* Take a step downhill.
1242 * In theory, we should minimize the function along this direction.
1243 * That is quite possible, but it turns out to take 5-10 function evaluations
1244 * for each line. However, we dont really need to find the exact minimum -
1245 * it is much better to start a new CG step in a modified direction as soon
1246 * as we are close to it. This will save a lot of energy evaluations.
1248 * In practice, we just try to take a single step.
1249 * If it worked (i.e. lowered the energy), we increase the stepsize but
1250 * the continue straight to the next CG step without trying to find any minimum.
1251 * If it didn't work (higher energy), there must be a minimum somewhere between
1252 * the old position and the new one.
1254 * Due to the finite numerical accuracy, it turns out that it is a good idea
1255 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1256 * This leads to lower final energies in the tests I've done. / Erik
1258 s_a->epot = s_min->epot;
1259 a = 0.0;
1260 c = a + stepsize; /* reference position along line is zero */
1262 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1264 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1265 s_min, top, mdAtoms, fr, vsite, constr,
1266 nrnb, wcycle);
1269 /* Take a trial step (new coords in s_c) */
1270 do_em_step(cr, inputrec, mdatoms, s_min, c, &s_min->s.cg_p, s_c,
1271 constr, -1);
1273 neval++;
1274 /* Calculate energy for the trial step */
1275 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1277 /* Calc derivative along line */
1278 const rvec *pc = as_rvec_array(s_c->s.cg_p.data());
1279 const rvec *sfc = as_rvec_array(s_c->f.data());
1280 double gpc = 0;
1281 for (int i = 0; i < mdatoms->homenr; i++)
1283 for (m = 0; m < DIM; m++)
1285 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1288 /* Sum the gradient along the line across CPUs */
1289 if (PAR(cr))
1291 gmx_sumd(1, &gpc, cr);
1294 /* This is the max amount of increase in energy we tolerate */
1295 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1297 /* Accept the step if the energy is lower, or if it is not significantly higher
1298 * and the line derivative is still negative.
1300 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1302 foundlower = TRUE;
1303 /* Great, we found a better energy. Increase step for next iteration
1304 * if we are still going down, decrease it otherwise
1306 if (gpc < 0)
1308 stepsize *= 1.618034; /* The golden section */
1310 else
1312 stepsize *= 0.618034; /* 1/golden section */
1315 else
1317 /* New energy is the same or higher. We will have to do some work
1318 * to find a smaller value in the interval. Take smaller step next time!
1320 foundlower = FALSE;
1321 stepsize *= 0.618034;
1327 /* OK, if we didn't find a lower value we will have to locate one now - there must
1328 * be one in the interval [a=0,c].
1329 * The same thing is valid here, though: Don't spend dozens of iterations to find
1330 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1331 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1333 * I also have a safeguard for potentially really pathological functions so we never
1334 * take more than 20 steps before we give up ...
1336 * If we already found a lower value we just skip this step and continue to the update.
1338 double gpb;
1339 if (!foundlower)
1341 nminstep = 0;
1345 /* Select a new trial point.
1346 * If the derivatives at points a & c have different sign we interpolate to zero,
1347 * otherwise just do a bisection.
1349 if (gpa < 0 && gpc > 0)
1351 b = a + gpa*(a-c)/(gpc-gpa);
1353 else
1355 b = 0.5*(a+c);
1358 /* safeguard if interpolation close to machine accuracy causes errors:
1359 * never go outside the interval
1361 if (b <= a || b >= c)
1363 b = 0.5*(a+c);
1366 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1368 /* Reload the old state */
1369 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1370 s_min, top, mdAtoms, fr, vsite, constr,
1371 nrnb, wcycle);
1374 /* Take a trial step to this new point - new coords in s_b */
1375 do_em_step(cr, inputrec, mdatoms, s_min, b, &s_min->s.cg_p, s_b,
1376 constr, -1);
1378 neval++;
1379 /* Calculate energy for the trial step */
1380 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1382 /* p does not change within a step, but since the domain decomposition
1383 * might change, we have to use cg_p of s_b here.
1385 const rvec *pb = as_rvec_array(s_b->s.cg_p.data());
1386 const rvec *sfb = as_rvec_array(s_b->f.data());
1387 gpb = 0;
1388 for (int i = 0; i < mdatoms->homenr; i++)
1390 for (m = 0; m < DIM; m++)
1392 gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1395 /* Sum the gradient along the line across CPUs */
1396 if (PAR(cr))
1398 gmx_sumd(1, &gpb, cr);
1401 if (debug)
1403 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1404 s_a->epot, s_b->epot, s_c->epot, gpb);
1407 epot_repl = s_b->epot;
1409 /* Keep one of the intervals based on the value of the derivative at the new point */
1410 if (gpb > 0)
1412 /* Replace c endpoint with b */
1413 swap_em_state(&s_b, &s_c);
1414 c = b;
1415 gpc = gpb;
1417 else
1419 /* Replace a endpoint with b */
1420 swap_em_state(&s_b, &s_a);
1421 a = b;
1422 gpa = gpb;
1426 * Stop search as soon as we find a value smaller than the endpoints.
1427 * Never run more than 20 steps, no matter what.
1429 nminstep++;
1431 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1432 (nminstep < 20));
1434 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1435 nminstep >= 20)
1437 /* OK. We couldn't find a significantly lower energy.
1438 * If beta==0 this was steepest descent, and then we give up.
1439 * If not, set beta=0 and restart with steepest descent before quitting.
1441 if (beta == 0.0)
1443 /* Converged */
1444 converged = TRUE;
1445 break;
1447 else
1449 /* Reset memory before giving up */
1450 beta = 0.0;
1451 continue;
1455 /* Select min energy state of A & C, put the best in B.
1457 if (s_c->epot < s_a->epot)
1459 if (debug)
1461 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1462 s_c->epot, s_a->epot);
1464 swap_em_state(&s_b, &s_c);
1465 gpb = gpc;
1467 else
1469 if (debug)
1471 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1472 s_a->epot, s_c->epot);
1474 swap_em_state(&s_b, &s_a);
1475 gpb = gpa;
1479 else
1481 if (debug)
1483 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1484 s_c->epot);
1486 swap_em_state(&s_b, &s_c);
1487 gpb = gpc;
1490 /* new search direction */
1491 /* beta = 0 means forget all memory and restart with steepest descents. */
1492 if (nstcg && ((step % nstcg) == 0))
1494 beta = 0.0;
1496 else
1498 /* s_min->fnorm cannot be zero, because then we would have converged
1499 * and broken out.
1502 /* Polak-Ribiere update.
1503 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1505 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1507 /* Limit beta to prevent oscillations */
1508 if (fabs(beta) > 5.0)
1510 beta = 0.0;
1514 /* update positions */
1515 swap_em_state(&s_min, &s_b);
1516 gpa = gpb;
1518 /* Print it if necessary */
1519 if (MASTER(cr))
1521 if (mdrunOptions.verbose)
1523 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1524 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1525 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1526 s_min->fmax, s_min->a_fmax+1);
1527 fflush(stderr);
1529 /* Store the new (lower) energies */
1530 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1531 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1532 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1534 do_log = do_per_step(step, inputrec->nstlog);
1535 do_ene = do_per_step(step, inputrec->nstenergy);
1537 /* Prepare IMD energy record, if bIMD is TRUE. */
1538 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1540 if (do_log)
1542 print_ebin_header(fplog, step, step);
1544 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1545 do_log ? fplog : nullptr, step, step, eprNORMAL,
1546 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1549 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1550 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
1552 IMD_send_positions(inputrec->imd);
1555 /* Stop when the maximum force lies below tolerance.
1556 * If we have reached machine precision, converged is already set to true.
1558 converged = converged || (s_min->fmax < inputrec->em_tol);
1560 } /* End of the loop */
1562 /* IMD cleanup, if bIMD is TRUE. */
1563 IMD_finalize(inputrec->bIMD, inputrec->imd);
1565 if (converged)
1567 step--; /* we never took that last step in this case */
1570 if (s_min->fmax > inputrec->em_tol)
1572 if (MASTER(cr))
1574 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1575 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1577 converged = FALSE;
1580 if (MASTER(cr))
1582 /* If we printed energy and/or logfile last step (which was the last step)
1583 * we don't have to do it again, but otherwise print the final values.
1585 if (!do_log)
1587 /* Write final value to log since we didn't do anything the last step */
1588 print_ebin_header(fplog, step, step);
1590 if (!do_ene || !do_log)
1592 /* Write final energy file entries */
1593 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1594 !do_log ? fplog : nullptr, step, step, eprNORMAL,
1595 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1599 /* Print some stuff... */
1600 if (MASTER(cr))
1602 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1605 /* IMPORTANT!
1606 * For accurate normal mode calculation it is imperative that we
1607 * store the last conformation into the full precision binary trajectory.
1609 * However, we should only do it if we did NOT already write this step
1610 * above (which we did if do_x or do_f was true).
1612 do_x = !do_per_step(step, inputrec->nstxout);
1613 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1615 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1616 top_global, inputrec, step,
1617 s_min, state_global, observablesHistory);
1620 if (MASTER(cr))
1622 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1623 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1624 s_min, sqrtNumAtoms);
1625 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1626 s_min, sqrtNumAtoms);
1628 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1631 finish_em(cr, outf, walltime_accounting, wcycle);
1633 /* To print the actual number of steps we needed somewhere */
1634 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1638 void
1639 Integrator::do_lbfgs()
1641 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1642 em_state_t ems;
1643 gmx_localtop_t *top;
1644 gmx_enerdata_t *enerd;
1645 gmx_global_stat_t gstat;
1646 t_graph *graph;
1647 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1648 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1649 real *rho, *alpha, *p, *s, **dx, **dg;
1650 real a, b, c, maxdelta, delta;
1651 real diag, Epot0;
1652 real dgdx, dgdg, sq, yr, beta;
1653 t_mdebin *mdebin;
1654 gmx_bool converged;
1655 rvec mu_tot;
1656 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1657 tensor vir, pres;
1658 int start, end, number_steps;
1659 gmx_mdoutf_t outf;
1660 int i, k, m, n, gf, step;
1661 int mdof_flags;
1662 auto mdatoms = mdAtoms->mdatoms();
1664 if (PAR(cr))
1666 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1669 if (nullptr != constr)
1671 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).");
1674 n = 3*state_global->natoms;
1675 nmaxcorr = inputrec->nbfgscorr;
1677 snew(frozen, n);
1679 snew(p, n);
1680 snew(rho, nmaxcorr);
1681 snew(alpha, nmaxcorr);
1683 snew(dx, nmaxcorr);
1684 for (i = 0; i < nmaxcorr; i++)
1686 snew(dx[i], n);
1689 snew(dg, nmaxcorr);
1690 for (i = 0; i < nmaxcorr; i++)
1692 snew(dg[i], n);
1695 step = 0;
1696 neval = 0;
1698 /* Init em */
1699 init_em(fplog, LBFGS, cr, ms, outputProvider, inputrec, mdrunOptions,
1700 state_global, top_global, &ems, &top,
1701 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1702 vsite, constr, nullptr,
1703 nfile, fnm, &outf, &mdebin, wcycle);
1705 start = 0;
1706 end = mdatoms->homenr;
1708 /* We need 4 working states */
1709 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1710 em_state_t *sa = &s0;
1711 em_state_t *sb = &s1;
1712 em_state_t *sc = &s2;
1713 em_state_t *last = &s3;
1714 /* Initialize by copying the state from ems (we could skip x and f here) */
1715 *sa = ems;
1716 *sb = ems;
1717 *sc = ems;
1719 /* Print to log file */
1720 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1722 do_log = do_ene = do_x = do_f = TRUE;
1724 /* Max number of steps */
1725 number_steps = inputrec->nsteps;
1727 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1728 gf = 0;
1729 for (i = start; i < end; i++)
1731 if (mdatoms->cFREEZE)
1733 gf = mdatoms->cFREEZE[i];
1735 for (m = 0; m < DIM; m++)
1737 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1740 if (MASTER(cr))
1742 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1744 if (fplog)
1746 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1749 if (vsite)
1751 construct_vsites(vsite, as_rvec_array(state_global->x.data()), 1, nullptr,
1752 top->idef.iparams, top->idef.il,
1753 fr->ePBC, fr->bMolPBC, cr, 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
1760 neval++;
1761 EnergyEvaluator energyEvaluator {
1762 fplog, cr, ms,
1763 top_global, top,
1764 inputrec, nrnb, wcycle, gstat,
1765 vsite, constr, fcd, graph,
1766 mdAtoms, fr, enerd
1768 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1770 if (MASTER(cr))
1772 /* Copy stuff to the energy bin for easy printing etc. */
1773 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1774 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
1775 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1777 print_ebin_header(fplog, step, step);
1778 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1779 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1782 /* Set the initial step.
1783 * since it will be multiplied by the non-normalized search direction
1784 * vector (force vector the first time), we scale it by the
1785 * norm of the force.
1788 if (MASTER(cr))
1790 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1791 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1792 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1793 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1794 fprintf(stderr, "\n");
1795 /* and copy to the log file too... */
1796 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1797 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1798 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1799 fprintf(fplog, "\n");
1802 // Point is an index to the memory of search directions, where 0 is the first one.
1803 point = 0;
1805 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1806 real *fInit = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1807 for (i = 0; i < n; i++)
1809 if (!frozen[i])
1811 dx[point][i] = fInit[i]; /* Initial search direction */
1813 else
1815 dx[point][i] = 0;
1819 // Stepsize will be modified during the search, and actually it is not critical
1820 // (the main efficiency in the algorithm comes from changing directions), but
1821 // we still need an initial value, so estimate it as the inverse of the norm
1822 // so we take small steps where the potential fluctuates a lot.
1823 stepsize = 1.0/ems.fnorm;
1825 /* Start the loop over BFGS steps.
1826 * Each successful step is counted, and we continue until
1827 * we either converge or reach the max number of steps.
1830 ncorr = 0;
1832 /* Set the gradient from the force */
1833 converged = FALSE;
1834 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1837 /* Write coordinates if necessary */
1838 do_x = do_per_step(step, inputrec->nstxout);
1839 do_f = do_per_step(step, inputrec->nstfout);
1841 mdof_flags = 0;
1842 if (do_x)
1844 mdof_flags |= MDOF_X;
1847 if (do_f)
1849 mdof_flags |= MDOF_F;
1852 if (inputrec->bIMD)
1854 mdof_flags |= MDOF_IMD;
1857 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1858 top_global, step, (real)step, &ems.s, state_global, 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 */
1863 s = dx[point];
1865 real *xx = static_cast<real *>(as_rvec_array(ems.s.x.data())[0]);
1866 real *ff = static_cast<real *>(as_rvec_array(ems.f.data())[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++)
1879 tmp = fabs(xx[i]);
1880 if (tmp < 1.0)
1882 tmp = 1.0;
1884 tmp = s[i]/tmp;
1885 minstep += tmp*tmp;
1887 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1889 if (stepsize < minstep)
1891 converged = TRUE;
1892 break;
1895 // Before taking any steps along the line, store the old position
1896 *last = ems;
1897 real *lastx = static_cast<real *>(as_rvec_array(last->s.x.data())[0]);
1898 real *lastf = static_cast<real *>(as_rvec_array(last->f.data())[0]);
1899 Epot0 = ems.epot;
1901 *sa = ems;
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
1930 a = 0.0;
1932 // Check stepsize first. We do not allow displacements
1933 // larger than emstep.
1937 // Pick a new position C by adding stepsize to A.
1938 c = a + stepsize;
1940 // Calculate what the largest change in any individual coordinate
1941 // would be (translation along line * gradient along line)
1942 maxdelta = 0;
1943 for (i = 0; i < n; i++)
1945 delta = c*s[i];
1946 if (delta > maxdelta)
1948 maxdelta = delta;
1951 // If any displacement is larger than the stepsize limit, reduce the step
1952 if (maxdelta > inputrec->em_stepsize)
1954 stepsize *= 0.1;
1957 while (maxdelta > inputrec->em_stepsize);
1959 // Take a trial step and move the coordinate array xc[] to position C
1960 real *xc = static_cast<real *>(as_rvec_array(sc->s.x.data())[0]);
1961 for (i = 0; i < n; i++)
1963 xc[i] = lastx[i] + c*s[i];
1966 neval++;
1967 // Calculate energy for the trial step in position C
1968 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
1970 // Calc line gradient in position C
1971 real *fc = static_cast<real *>(as_rvec_array(sc->f.data())[0]);
1972 for (gpc = 0, i = 0; i < n; i++)
1974 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1976 /* Sum the gradient along the line across CPUs */
1977 if (PAR(cr))
1979 gmx_sumd(1, &gpc, cr);
1982 // This is the max amount of increase in energy we tolerate.
1983 // By allowing VERY small changes (close to numerical precision) we
1984 // frequently find even better (lower) final energies.
1985 tmp = sqrt(GMX_REAL_EPS)*fabs(sa->epot);
1987 // Accept the step if the energy is lower in the new position C (compared to A),
1988 // or if it is not significantly higher and the line derivative is still negative.
1989 if (sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp)))
1991 // Great, we found a better energy. We no longer try to alter the
1992 // stepsize, but simply accept this new better position. The we select a new
1993 // search direction instead, which will be much more efficient than continuing
1994 // to take smaller steps along a line. Set fnorm based on the new C position,
1995 // which will be used to update the stepsize to 1/fnorm further down.
1996 foundlower = TRUE;
1998 else
2000 // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2001 // or higher than in point A. In this case it is pointless to move to point C,
2002 // so we will have to do more iterations along the same line to find a smaller
2003 // value in the interval [A=0.0,C].
2004 // Here, A is still 0.0, but that will change when we do a search in the interval
2005 // [0.0,C] below. That search we will do by interpolation or bisection rather
2006 // than with the stepsize, so no need to modify it. For the next search direction
2007 // it will be reset to 1/fnorm anyway.
2008 foundlower = FALSE;
2011 if (!foundlower)
2013 // OK, if we didn't find a lower value we will have to locate one now - there must
2014 // be one in the interval [a,c].
2015 // The same thing is valid here, though: Don't spend dozens of iterations to find
2016 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2017 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2018 // I also have a safeguard for potentially really pathological functions so we never
2019 // take more than 20 steps before we give up.
2020 // If we already found a lower value we just skip this step and continue to the update.
2021 real fnorm = 0;
2022 nminstep = 0;
2025 // Select a new trial point B in the interval [A,C].
2026 // If the derivatives at points a & c have different sign we interpolate to zero,
2027 // otherwise just do a bisection since there might be multiple minima/maxima
2028 // inside the interval.
2029 if (gpa < 0 && gpc > 0)
2031 b = a + gpa*(a-c)/(gpc-gpa);
2033 else
2035 b = 0.5*(a+c);
2038 /* safeguard if interpolation close to machine accuracy causes errors:
2039 * never go outside the interval
2041 if (b <= a || b >= c)
2043 b = 0.5*(a+c);
2046 // Take a trial step to point B
2047 real *xb = static_cast<real *>(as_rvec_array(sb->s.x.data())[0]);
2048 for (i = 0; i < n; i++)
2050 xb[i] = lastx[i] + b*s[i];
2053 neval++;
2054 // Calculate energy for the trial step in point B
2055 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2056 fnorm = sb->fnorm;
2058 // Calculate gradient in point B
2059 real *fb = static_cast<real *>(as_rvec_array(sb->f.data())[0]);
2060 for (gpb = 0, i = 0; i < n; i++)
2062 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2065 /* Sum the gradient along the line across CPUs */
2066 if (PAR(cr))
2068 gmx_sumd(1, &gpb, cr);
2071 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2072 // at the new point B, and rename the endpoints of this new interval A and C.
2073 if (gpb > 0)
2075 /* Replace c endpoint with b */
2076 c = b;
2077 /* swap states b and c */
2078 swap_em_state(&sb, &sc);
2080 else
2082 /* Replace a endpoint with b */
2083 a = b;
2084 /* swap states a and b */
2085 swap_em_state(&sa, &sb);
2089 * Stop search as soon as we find a value smaller than the endpoints,
2090 * or if the tolerance is below machine precision.
2091 * Never run more than 20 steps, no matter what.
2093 nminstep++;
2095 while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2097 if (fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2099 /* OK. We couldn't find a significantly lower energy.
2100 * If ncorr==0 this was steepest descent, and then we give up.
2101 * If not, reset memory to restart as steepest descent before quitting.
2103 if (ncorr == 0)
2105 /* Converged */
2106 converged = TRUE;
2107 break;
2109 else
2111 /* Reset memory */
2112 ncorr = 0;
2113 /* Search in gradient direction */
2114 for (i = 0; i < n; i++)
2116 dx[point][i] = ff[i];
2118 /* Reset stepsize */
2119 stepsize = 1.0/fnorm;
2120 continue;
2124 /* Select min energy state of A & C, put the best in xx/ff/Epot
2126 if (sc->epot < sa->epot)
2128 /* Use state C */
2129 ems = *sc;
2130 step_taken = c;
2132 else
2134 /* Use state A */
2135 ems = *sa;
2136 step_taken = a;
2140 else
2142 /* found lower */
2143 /* Use state C */
2144 ems = *sc;
2145 step_taken = c;
2148 /* Update the memory information, and calculate a new
2149 * approximation of the inverse hessian
2152 /* Have new data in Epot, xx, ff */
2153 if (ncorr < nmaxcorr)
2155 ncorr++;
2158 for (i = 0; i < n; i++)
2160 dg[point][i] = lastf[i]-ff[i];
2161 dx[point][i] *= step_taken;
2164 dgdg = 0;
2165 dgdx = 0;
2166 for (i = 0; i < n; i++)
2168 dgdg += dg[point][i]*dg[point][i];
2169 dgdx += dg[point][i]*dx[point][i];
2172 diag = dgdx/dgdg;
2174 rho[point] = 1.0/dgdx;
2175 point++;
2177 if (point >= nmaxcorr)
2179 point = 0;
2182 /* Update */
2183 for (i = 0; i < n; i++)
2185 p[i] = ff[i];
2188 cp = point;
2190 /* Recursive update. First go back over the memory points */
2191 for (k = 0; k < ncorr; k++)
2193 cp--;
2194 if (cp < 0)
2196 cp = ncorr-1;
2199 sq = 0;
2200 for (i = 0; i < n; i++)
2202 sq += dx[cp][i]*p[i];
2205 alpha[cp] = rho[cp]*sq;
2207 for (i = 0; i < n; i++)
2209 p[i] -= alpha[cp]*dg[cp][i];
2213 for (i = 0; i < n; i++)
2215 p[i] *= diag;
2218 /* And then go forward again */
2219 for (k = 0; k < ncorr; k++)
2221 yr = 0;
2222 for (i = 0; i < n; i++)
2224 yr += p[i]*dg[cp][i];
2227 beta = rho[cp]*yr;
2228 beta = alpha[cp]-beta;
2230 for (i = 0; i < n; i++)
2232 p[i] += beta*dx[cp][i];
2235 cp++;
2236 if (cp >= ncorr)
2238 cp = 0;
2242 for (i = 0; i < n; i++)
2244 if (!frozen[i])
2246 dx[point][i] = p[i];
2248 else
2250 dx[point][i] = 0;
2254 /* Print it if necessary */
2255 if (MASTER(cr))
2257 if (mdrunOptions.verbose)
2259 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2260 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2261 step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2262 fflush(stderr);
2264 /* Store the new (lower) energies */
2265 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2266 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
2267 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2268 do_log = do_per_step(step, inputrec->nstlog);
2269 do_ene = do_per_step(step, inputrec->nstenergy);
2270 if (do_log)
2272 print_ebin_header(fplog, step, step);
2274 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2275 do_log ? fplog : nullptr, step, step, eprNORMAL,
2276 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2279 /* Send x and E to IMD client, if bIMD is TRUE. */
2280 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
2282 IMD_send_positions(inputrec->imd);
2285 // Reset stepsize in we are doing more iterations
2286 stepsize = 1.0/ems.fnorm;
2288 /* Stop when the maximum force lies below tolerance.
2289 * If we have reached machine precision, converged is already set to true.
2291 converged = converged || (ems.fmax < inputrec->em_tol);
2293 } /* End of the loop */
2295 /* IMD cleanup, if bIMD is TRUE. */
2296 IMD_finalize(inputrec->bIMD, inputrec->imd);
2298 if (converged)
2300 step--; /* we never took that last step in this case */
2303 if (ems.fmax > inputrec->em_tol)
2305 if (MASTER(cr))
2307 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2308 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2310 converged = FALSE;
2313 /* If we printed energy and/or logfile last step (which was the last step)
2314 * we don't have to do it again, but otherwise print the final values.
2316 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2318 print_ebin_header(fplog, step, step);
2320 if (!do_ene || !do_log) /* Write final energy file entries */
2322 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2323 !do_log ? fplog : nullptr, step, step, eprNORMAL,
2324 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2327 /* Print some stuff... */
2328 if (MASTER(cr))
2330 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2333 /* IMPORTANT!
2334 * For accurate normal mode calculation it is imperative that we
2335 * store the last conformation into the full precision binary trajectory.
2337 * However, we should only do it if we did NOT already write this step
2338 * above (which we did if do_x or do_f was true).
2340 do_x = !do_per_step(step, inputrec->nstxout);
2341 do_f = !do_per_step(step, inputrec->nstfout);
2342 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2343 top_global, inputrec, step,
2344 &ems, state_global, observablesHistory);
2346 if (MASTER(cr))
2348 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2349 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2350 number_steps, &ems, sqrtNumAtoms);
2351 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2352 number_steps, &ems, sqrtNumAtoms);
2354 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2357 finish_em(cr, outf, walltime_accounting, wcycle);
2359 /* To print the actual number of steps we needed somewhere */
2360 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2363 void
2364 Integrator::do_steep()
2366 const char *SD = "Steepest Descents";
2367 gmx_localtop_t *top;
2368 gmx_enerdata_t *enerd;
2369 gmx_global_stat_t gstat;
2370 t_graph *graph;
2371 real stepsize;
2372 real ustep;
2373 gmx_mdoutf_t outf;
2374 t_mdebin *mdebin;
2375 gmx_bool bDone, bAbort, do_x, do_f;
2376 tensor vir, pres;
2377 rvec mu_tot;
2378 int nsteps;
2379 int count = 0;
2380 int steps_accepted = 0;
2381 auto mdatoms = mdAtoms->mdatoms();
2383 /* Create 2 states on the stack and extract pointers that we will swap */
2384 em_state_t s0 {}, s1 {};
2385 em_state_t *s_min = &s0;
2386 em_state_t *s_try = &s1;
2388 /* Init em and store the local state in s_try */
2389 init_em(fplog, SD, cr, ms, outputProvider, inputrec, mdrunOptions,
2390 state_global, top_global, s_try, &top,
2391 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2392 vsite, constr, nullptr,
2393 nfile, fnm, &outf, &mdebin, wcycle);
2395 /* Print to log file */
2396 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2398 /* Set variables for stepsize (in nm). This is the largest
2399 * step that we are going to make in any direction.
2401 ustep = inputrec->em_stepsize;
2402 stepsize = 0;
2404 /* Max number of steps */
2405 nsteps = inputrec->nsteps;
2407 if (MASTER(cr))
2409 /* Print to the screen */
2410 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2412 if (fplog)
2414 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2416 EnergyEvaluator energyEvaluator {
2417 fplog, cr, ms,
2418 top_global, top,
2419 inputrec, nrnb, wcycle, gstat,
2420 vsite, constr, fcd, graph,
2421 mdAtoms, fr, enerd
2424 /**** HERE STARTS THE LOOP ****
2425 * count is the counter for the number of steps
2426 * bDone will be TRUE when the minimization has converged
2427 * bAbort will be TRUE when nsteps steps have been performed or when
2428 * the stepsize becomes smaller than is reasonable for machine precision
2430 count = 0;
2431 bDone = FALSE;
2432 bAbort = FALSE;
2433 while (!bDone && !bAbort)
2435 bAbort = (nsteps >= 0) && (count == nsteps);
2437 /* set new coordinates, except for first step */
2438 bool validStep = true;
2439 if (count > 0)
2441 validStep =
2442 do_em_step(cr, inputrec, mdatoms,
2443 s_min, stepsize, &s_min->f, s_try,
2444 constr, count);
2447 if (validStep)
2449 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2451 else
2453 // Signal constraint error during stepping with energy=inf
2454 s_try->epot = std::numeric_limits<real>::infinity();
2457 if (MASTER(cr))
2459 print_ebin_header(fplog, count, count);
2462 if (count == 0)
2464 s_min->epot = s_try->epot;
2467 /* Print it if necessary */
2468 if (MASTER(cr))
2470 if (mdrunOptions.verbose)
2472 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2473 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2474 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2475 fflush(stderr);
2478 if ( (count == 0) || (s_try->epot < s_min->epot) )
2480 /* Store the new (lower) energies */
2481 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2482 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2483 s_try->s.box, nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2485 /* Prepare IMD energy record, if bIMD is TRUE. */
2486 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2488 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2489 do_per_step(steps_accepted, inputrec->nstdisreout),
2490 do_per_step(steps_accepted, inputrec->nstorireout),
2491 fplog, count, count, eprNORMAL,
2492 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2493 fflush(fplog);
2497 /* Now if the new energy is smaller than the previous...
2498 * or if this is the first step!
2499 * or if we did random steps!
2502 if ( (count == 0) || (s_try->epot < s_min->epot) )
2504 steps_accepted++;
2506 /* Test whether the convergence criterion is met... */
2507 bDone = (s_try->fmax < inputrec->em_tol);
2509 /* Copy the arrays for force, positions and energy */
2510 /* The 'Min' array always holds the coords and forces of the minimal
2511 sampled energy */
2512 swap_em_state(&s_min, &s_try);
2513 if (count > 0)
2515 ustep *= 1.2;
2518 /* Write to trn, if necessary */
2519 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2520 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2521 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2522 top_global, inputrec, count,
2523 s_min, state_global, observablesHistory);
2525 else
2527 /* If energy is not smaller make the step smaller... */
2528 ustep *= 0.5;
2530 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2532 /* Reload the old state */
2533 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2534 s_min, top, mdAtoms, fr, vsite, constr,
2535 nrnb, wcycle);
2539 /* Determine new step */
2540 stepsize = ustep/s_min->fmax;
2542 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2543 #if GMX_DOUBLE
2544 if (count == nsteps || ustep < 1e-12)
2545 #else
2546 if (count == nsteps || ustep < 1e-6)
2547 #endif
2549 if (MASTER(cr))
2551 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != nullptr);
2552 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != nullptr);
2554 bAbort = TRUE;
2557 /* Send IMD energies and positions, if bIMD is TRUE. */
2558 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box,
2559 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
2560 inputrec, 0, wcycle) &&
2561 MASTER(cr))
2563 IMD_send_positions(inputrec->imd);
2566 count++;
2567 } /* End of the loop */
2569 /* IMD cleanup, if bIMD is TRUE. */
2570 IMD_finalize(inputrec->bIMD, inputrec->imd);
2572 /* Print some data... */
2573 if (MASTER(cr))
2575 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2577 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2578 top_global, inputrec, count,
2579 s_min, state_global, observablesHistory);
2581 if (MASTER(cr))
2583 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2585 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2586 s_min, sqrtNumAtoms);
2587 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2588 s_min, sqrtNumAtoms);
2591 finish_em(cr, outf, walltime_accounting, wcycle);
2593 /* To print the actual number of steps we needed somewhere */
2594 inputrec->nsteps = count;
2596 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2599 void
2600 Integrator::do_nm()
2602 const char *NM = "Normal Mode Analysis";
2603 gmx_mdoutf_t outf;
2604 int nnodes, node;
2605 gmx_localtop_t *top;
2606 gmx_enerdata_t *enerd;
2607 gmx_global_stat_t gstat;
2608 t_graph *graph;
2609 tensor vir, pres;
2610 rvec mu_tot;
2611 rvec *fneg, *dfdx;
2612 gmx_bool bSparse; /* use sparse matrix storage format */
2613 size_t sz;
2614 gmx_sparsematrix_t * sparse_matrix = nullptr;
2615 real * full_matrix = nullptr;
2617 /* added with respect to mdrun */
2618 int row, col;
2619 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2620 real x_min;
2621 bool bIsMaster = MASTER(cr);
2622 auto mdatoms = mdAtoms->mdatoms();
2624 if (constr != nullptr)
2626 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2629 gmx_shellfc_t *shellfc;
2631 em_state_t state_work {};
2633 /* Init em and store the local state in state_minimum */
2634 init_em(fplog, NM, cr, ms, outputProvider, inputrec, mdrunOptions,
2635 state_global, top_global, &state_work, &top,
2636 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2637 vsite, constr, &shellfc,
2638 nfile, fnm, &outf, nullptr, wcycle);
2640 std::vector<size_t> atom_index = get_atom_index(top_global);
2641 snew(fneg, atom_index.size());
2642 snew(dfdx, atom_index.size());
2644 #if !GMX_DOUBLE
2645 if (bIsMaster)
2647 fprintf(stderr,
2648 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2649 " which MIGHT not be accurate enough for normal mode analysis.\n"
2650 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2651 " are fairly modest even if you recompile in double precision.\n\n");
2653 #endif
2655 /* Check if we can/should use sparse storage format.
2657 * Sparse format is only useful when the Hessian itself is sparse, which it
2658 * will be when we use a cutoff.
2659 * For small systems (n<1000) it is easier to always use full matrix format, though.
2661 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2663 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2664 bSparse = FALSE;
2666 else if (atom_index.size() < 1000)
2668 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%d), using full Hessian format.",
2669 atom_index.size());
2670 bSparse = FALSE;
2672 else
2674 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2675 bSparse = TRUE;
2678 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2679 sz = DIM*atom_index.size();
2681 fprintf(stderr, "Allocating Hessian memory...\n\n");
2683 if (bSparse)
2685 sparse_matrix = gmx_sparsematrix_init(sz);
2686 sparse_matrix->compressed_symmetric = TRUE;
2688 else
2690 snew(full_matrix, sz*sz);
2693 init_nrnb(nrnb);
2696 /* Write start time and temperature */
2697 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2699 /* fudge nr of steps to nr of atoms */
2700 inputrec->nsteps = atom_index.size()*2;
2702 if (bIsMaster)
2704 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2705 *(top_global->name), (int)inputrec->nsteps);
2708 nnodes = cr->nnodes;
2710 /* Make evaluate_energy do a single node force calculation */
2711 cr->nnodes = 1;
2712 EnergyEvaluator energyEvaluator {
2713 fplog, cr, ms,
2714 top_global, top,
2715 inputrec, nrnb, wcycle, gstat,
2716 vsite, constr, fcd, graph,
2717 mdAtoms, fr, enerd
2719 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2720 cr->nnodes = nnodes;
2722 /* if forces are not small, warn user */
2723 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2725 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2726 if (state_work.fmax > 1.0e-3)
2728 GMX_LOG(mdlog.warning).appendText(
2729 "The force is probably not small enough to "
2730 "ensure that you are at a minimum.\n"
2731 "Be aware that negative eigenvalues may occur\n"
2732 "when the resulting matrix is diagonalized.");
2735 /***********************************************************
2737 * Loop over all pairs in matrix
2739 * do_force called twice. Once with positive and
2740 * once with negative displacement
2742 ************************************************************/
2744 /* Steps are divided one by one over the nodes */
2745 bool bNS = true;
2746 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2748 size_t atom = atom_index[aid];
2749 for (size_t d = 0; d < DIM; d++)
2751 gmx_int64_t step = 0;
2752 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2753 double t = 0;
2755 x_min = state_work.s.x[atom][d];
2757 for (unsigned int dx = 0; (dx < 2); dx++)
2759 if (dx == 0)
2761 state_work.s.x[atom][d] = x_min - der_range;
2763 else
2765 state_work.s.x[atom][d] = x_min + der_range;
2768 /* Make evaluate_energy do a single node force calculation */
2769 cr->nnodes = 1;
2770 if (shellfc)
2772 /* Now is the time to relax the shells */
2773 relax_shell_flexcon(fplog,
2776 mdrunOptions.verbose,
2777 step,
2778 inputrec,
2779 bNS,
2780 force_flags,
2781 top,
2782 constr,
2783 enerd,
2784 fcd,
2785 &state_work.s,
2786 state_work.f,
2787 vir,
2788 mdatoms,
2789 nrnb,
2790 wcycle,
2791 graph,
2792 &top_global->groups,
2793 shellfc,
2796 mu_tot,
2797 vsite,
2798 DdOpenBalanceRegionBeforeForceComputation::no,
2799 DdCloseBalanceRegionAfterForceComputation::no);
2800 bNS = false;
2801 step++;
2803 else
2805 energyEvaluator.run(&state_work, mu_tot, vir, pres, atom*2+dx, FALSE);
2808 cr->nnodes = nnodes;
2810 if (dx == 0)
2812 for (size_t i = 0; i < atom_index.size(); i++)
2814 copy_rvec(state_work.f[atom_index[i]], fneg[i]);
2819 /* x is restored to original */
2820 state_work.s.x[atom][d] = x_min;
2822 for (size_t j = 0; j < atom_index.size(); j++)
2824 for (size_t k = 0; (k < DIM); k++)
2826 dfdx[j][k] =
2827 -(state_work.f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2831 if (!bIsMaster)
2833 #if GMX_MPI
2834 #define mpi_type GMX_MPI_REAL
2835 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2836 cr->nodeid, cr->mpi_comm_mygroup);
2837 #endif
2839 else
2841 for (node = 0; (node < nnodes && atom+node < atom_index.size()); node++)
2843 if (node > 0)
2845 #if GMX_MPI
2846 MPI_Status stat;
2847 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2848 cr->mpi_comm_mygroup, &stat);
2849 #undef mpi_type
2850 #endif
2853 row = (atom + node)*DIM + d;
2855 for (size_t j = 0; j < atom_index.size(); j++)
2857 for (size_t k = 0; k < DIM; k++)
2859 col = j*DIM + k;
2861 if (bSparse)
2863 if (col >= row && dfdx[j][k] != 0.0)
2865 gmx_sparsematrix_increment_value(sparse_matrix,
2866 row, col, dfdx[j][k]);
2869 else
2871 full_matrix[row*sz+col] = dfdx[j][k];
2878 if (mdrunOptions.verbose && fplog)
2880 fflush(fplog);
2883 /* write progress */
2884 if (bIsMaster && mdrunOptions.verbose)
2886 fprintf(stderr, "\rFinished step %d out of %d",
2887 static_cast<int>(std::min(atom+nnodes, atom_index.size())),
2888 static_cast<int>(atom_index.size()));
2889 fflush(stderr);
2893 if (bIsMaster)
2895 fprintf(stderr, "\n\nWriting Hessian...\n");
2896 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2899 finish_em(cr, outf, walltime_accounting, wcycle);
2901 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);
2904 } // namespace gmx