More const correctness
[gromacs.git] / src / gromacs / mdrun / minimize.cpp
blob33ee89e395c23216199d0e6425aa9507a1f0585b
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->gatindex[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 gmx::n_flexible_constraints(constr),
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 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 if (ir->eConstrAlg == econtSHAKE &&
429 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
431 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
432 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
435 if (!DOMAINDECOMP(cr))
437 set_constraints(constr, *top, ir, mdatoms, cr);
440 if (!ir->bContinuation)
442 /* Constrain the starting coordinates */
443 dvdl_constr = 0;
444 constrain(PAR(cr) ? nullptr : fplog, TRUE, TRUE, constr, &(*top)->idef,
445 ir, cr, ms, -1, 0, 1.0, mdatoms,
446 as_rvec_array(ems->s.x.data()),
447 as_rvec_array(ems->s.x.data()),
448 nullptr,
449 fr->bMolPBC, ems->s.box,
450 ems->s.lambda[efptFEP], &dvdl_constr,
451 nullptr, nullptr, nrnb, gmx::econqCoord);
455 if (PAR(cr))
457 *gstat = global_stat_init(ir);
459 else
461 *gstat = nullptr;
464 *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, nullptr, wcycle);
466 snew(*enerd, 1);
467 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
468 *enerd);
470 if (mdebin != nullptr)
472 /* Init bin for energy stuff */
473 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, nullptr);
476 clear_rvec(mu_tot);
477 calc_shifts(ems->s.box, fr->shift_vec);
480 //! Finalize the minimization
481 static void finish_em(const t_commrec *cr, gmx_mdoutf_t outf,
482 gmx_walltime_accounting_t walltime_accounting,
483 gmx_wallcycle_t wcycle)
485 if (!thisRankHasDuty(cr, DUTY_PME))
487 /* Tell the PME only node to finish */
488 gmx_pme_send_finish(cr);
491 done_mdoutf(outf);
493 em_time_end(walltime_accounting, wcycle);
496 //! Swap two different EM states during minimization
497 static void swap_em_state(em_state_t **ems1, em_state_t **ems2)
499 em_state_t *tmp;
501 tmp = *ems1;
502 *ems1 = *ems2;
503 *ems2 = tmp;
506 //! Save the EM trajectory
507 static void write_em_traj(FILE *fplog, const t_commrec *cr,
508 gmx_mdoutf_t outf,
509 gmx_bool bX, gmx_bool bF, const char *confout,
510 gmx_mtop_t *top_global,
511 t_inputrec *ir, gmx_int64_t step,
512 em_state_t *state,
513 t_state *state_global,
514 ObservablesHistory *observablesHistory)
516 int mdof_flags = 0;
518 if (bX)
520 mdof_flags |= MDOF_X;
522 if (bF)
524 mdof_flags |= MDOF_F;
527 /* If we want IMD output, set appropriate MDOF flag */
528 if (ir->bIMD)
530 mdof_flags |= MDOF_IMD;
533 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
534 top_global, step, (double)step,
535 &state->s, state_global, observablesHistory,
536 state->f);
538 if (confout != nullptr && MASTER(cr))
540 GMX_RELEASE_ASSERT(bX, "The code below assumes that (with domain decomposition), x is collected to state_global in the call above.");
541 /* With domain decomposition the call above collected the state->s.x
542 * into state_global->x. Without DD we copy the local state pointer.
544 if (!DOMAINDECOMP(cr))
546 state_global = &state->s;
549 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
551 /* Make molecules whole only for confout writing */
552 do_pbc_mtop(fplog, ir->ePBC, state->s.box, top_global,
553 as_rvec_array(state_global->x.data()));
556 write_sto_conf_mtop(confout,
557 *top_global->name, top_global,
558 as_rvec_array(state_global->x.data()), nullptr, ir->ePBC, state->s.box);
562 //! \brief Do one minimization step
564 // \returns true when the step succeeded, false when a constraint error occurred
565 static bool do_em_step(const t_commrec *cr,
566 const gmx_multisim_t *ms,
567 t_inputrec *ir, t_mdatoms *md,
568 gmx_bool bMolPBC,
569 em_state_t *ems1, real a, const PaddedRVecVector *force,
570 em_state_t *ems2,
571 gmx::Constraints *constr, gmx_localtop_t *top,
572 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
573 gmx_int64_t count)
576 t_state *s1, *s2;
577 int start, end;
578 real dvdl_constr;
579 int nthreads gmx_unused;
581 bool validStep = true;
583 s1 = &ems1->s;
584 s2 = &ems2->s;
586 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
588 gmx_incons("state mismatch in do_em_step");
591 s2->flags = s1->flags;
593 if (s2->natoms != s1->natoms)
595 state_change_natoms(s2, s1->natoms);
596 /* We need to allocate one element extra, since we might use
597 * (unaligned) 4-wide SIMD loads to access rvec entries.
599 ems2->f.resize(gmx::paddedRVecVectorSize(s2->natoms));
601 if (DOMAINDECOMP(cr) && s2->cg_gl.size() != s1->cg_gl.size())
603 s2->cg_gl.resize(s1->cg_gl.size());
606 copy_mat(s1->box, s2->box);
607 /* Copy free energy state */
608 s2->lambda = s1->lambda;
609 copy_mat(s1->box, s2->box);
611 start = 0;
612 end = md->homenr;
614 // cppcheck-suppress unreadVariable
615 nthreads = gmx_omp_nthreads_get(emntUpdate);
616 #pragma omp parallel num_threads(nthreads)
618 const rvec *x1 = as_rvec_array(s1->x.data());
619 rvec *x2 = as_rvec_array(s2->x.data());
620 const rvec *f = as_rvec_array(force->data());
622 int gf = 0;
623 #pragma omp for schedule(static) nowait
624 for (int i = start; i < end; i++)
628 if (md->cFREEZE)
630 gf = md->cFREEZE[i];
632 for (int m = 0; m < DIM; m++)
634 if (ir->opts.nFreeze[gf][m])
636 x2[i][m] = x1[i][m];
638 else
640 x2[i][m] = x1[i][m] + a*f[i][m];
644 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
647 if (s2->flags & (1<<estCGP))
649 /* Copy the CG p vector */
650 const rvec *p1 = as_rvec_array(s1->cg_p.data());
651 rvec *p2 = as_rvec_array(s2->cg_p.data());
652 #pragma omp for schedule(static) nowait
653 for (int i = start; i < end; i++)
655 // Trivial OpenMP block that does not throw
656 copy_rvec(p1[i], p2[i]);
660 if (DOMAINDECOMP(cr))
662 s2->ddp_count = s1->ddp_count;
664 /* OpenMP does not supported unsigned loop variables */
665 #pragma omp for schedule(static) nowait
666 for (int i = 0; i < static_cast<int>(s2->cg_gl.size()); i++)
668 s2->cg_gl[i] = s1->cg_gl[i];
670 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
674 if (constr)
676 wallcycle_start(wcycle, ewcCONSTR);
677 dvdl_constr = 0;
678 validStep =
679 constrain(nullptr, TRUE, TRUE, constr, &top->idef,
680 ir, cr, ms, count, 0, 1.0, md,
681 as_rvec_array(s1->x.data()), as_rvec_array(s2->x.data()),
682 nullptr, bMolPBC, s2->box,
683 s2->lambda[efptBONDED], &dvdl_constr,
684 nullptr, nullptr, nrnb, gmx::econqCoord);
685 wallcycle_stop(wcycle, ewcCONSTR);
687 // We should move this check to the different minimizers
688 if (!validStep && ir->eI != eiSteep)
690 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
691 EI(ir->eI), EI(eiSteep), EI(ir->eI));
695 return validStep;
698 //! Prepare EM for using domain decomposition parallellization
699 static void em_dd_partition_system(FILE *fplog, int step, const t_commrec *cr,
700 gmx_mtop_t *top_global, t_inputrec *ir,
701 em_state_t *ems, gmx_localtop_t *top,
702 gmx::MDAtoms *mdAtoms, t_forcerec *fr,
703 gmx_vsite_t *vsite, gmx::Constraints *constr,
704 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
706 /* Repartition the domain decomposition */
707 dd_partition_system(fplog, step, cr, FALSE, 1,
708 nullptr, top_global, ir,
709 &ems->s, &ems->f,
710 mdAtoms, top, fr, vsite, constr,
711 nrnb, wcycle, FALSE);
712 dd_store_state(cr->dd, &ems->s);
715 namespace
718 /*! \brief Class to handle the work of setting and doing an energy evaluation.
720 * This class is a mere aggregate of parameters to pass to evaluate an
721 * energy, so that future changes to names and types of them consume
722 * less time when refactoring other code.
724 * Aggregate initialization is used, for which the chief risk is that
725 * if a member is added at the end and not all initializer lists are
726 * updated, then the member will be value initialized, which will
727 * typically mean initialization to zero.
729 * We only want to construct one of these with an initializer list, so
730 * we explicitly delete the default constructor. */
731 class EnergyEvaluator
733 public:
734 //! We only intend to construct such objects with an initializer list.
735 #if __GNUC__ > 4 || (__GNUC__ == 4 && __GNUC_MINOR__ >= 9)
736 // Aspects of the C++11 spec changed after GCC 4.8.5, and
737 // compilation of the initializer list construction in
738 // runner.cpp fails in GCC 4.8.5.
739 EnergyEvaluator() = delete;
740 #endif
741 /*! \brief Evaluates an energy on the state in \c ems.
743 * \todo In practice, the same objects mu_tot, vir, and pres
744 * are always passed to this function, so we would rather have
745 * them as data members. However, their C-array types are
746 * unsuited for aggregate initialization. When the types
747 * improve, the call signature of this method can be reduced.
749 void run(em_state_t *ems, rvec mu_tot,
750 tensor vir, tensor pres,
751 gmx_int64_t count, gmx_bool bFirst);
752 //! Handles logging.
753 FILE *fplog;
754 //! Handles communication.
755 const t_commrec *cr;
756 //! Coordinates multi-simulations.
757 const gmx_multisim_t *ms;
758 //! Holds the simulation topology.
759 gmx_mtop_t *top_global;
760 //! Holds the domain topology.
761 gmx_localtop_t *top;
762 //! User input options.
763 t_inputrec *inputrec;
764 //! Manages flop accounting.
765 t_nrnb *nrnb;
766 //! Manages wall cycle accounting.
767 gmx_wallcycle_t wcycle;
768 //! Coordinates global reduction.
769 gmx_global_stat_t gstat;
770 //! Handles virtual sites.
771 gmx_vsite_t *vsite;
772 //! Handles constraints.
773 gmx::Constraints *constr;
774 //! Handles strange things.
775 t_fcdata *fcd;
776 //! Molecular graph for SHAKE.
777 t_graph *graph;
778 //! Per-atom data for this domain.
779 gmx::MDAtoms *mdAtoms;
780 //! Handles how to calculate the forces.
781 t_forcerec *fr;
782 //! Stores the computed energies.
783 gmx_enerdata_t *enerd;
786 void
787 EnergyEvaluator::run(em_state_t *ems, rvec mu_tot,
788 tensor vir, tensor pres,
789 gmx_int64_t count, gmx_bool bFirst)
791 real t;
792 gmx_bool bNS;
793 tensor force_vir, shake_vir, ekin;
794 real dvdl_constr, prescorr, enercorr, dvdlcorr;
795 real terminate = 0;
797 /* Set the time to the initial time, the time does not change during EM */
798 t = inputrec->init_t;
800 if (bFirst ||
801 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
803 /* This is the first state or an old state used before the last ns */
804 bNS = TRUE;
806 else
808 bNS = FALSE;
809 if (inputrec->nstlist > 0)
811 bNS = TRUE;
815 if (vsite)
817 construct_vsites(vsite, as_rvec_array(ems->s.x.data()), 1, nullptr,
818 top->idef.iparams, top->idef.il,
819 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
822 if (DOMAINDECOMP(cr) && bNS)
824 /* Repartition the domain decomposition */
825 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
826 ems, top, mdAtoms, fr, vsite, constr,
827 nrnb, wcycle);
830 /* Calc force & energy on new trial position */
831 /* do_force always puts the charge groups in the box and shifts again
832 * We do not unshift, so molecules are always whole in congrad.c
834 do_force(fplog, cr, ms, inputrec,
835 count, nrnb, wcycle, top, &top_global->groups,
836 ems->s.box, ems->s.x, &ems->s.hist,
837 ems->f, force_vir, mdAtoms->mdatoms(), enerd, fcd,
838 ems->s.lambda, graph, fr, vsite, mu_tot, t, nullptr,
839 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
840 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
841 (bNS ? GMX_FORCE_NS : 0),
842 DOMAINDECOMP(cr) ?
843 DdOpenBalanceRegionBeforeForceComputation::yes :
844 DdOpenBalanceRegionBeforeForceComputation::no,
845 DOMAINDECOMP(cr) ?
846 DdCloseBalanceRegionAfterForceComputation::yes :
847 DdCloseBalanceRegionAfterForceComputation::no);
849 /* Clear the unused shake virial and pressure */
850 clear_mat(shake_vir);
851 clear_mat(pres);
853 /* Communicate stuff when parallel */
854 if (PAR(cr) && inputrec->eI != eiNM)
856 wallcycle_start(wcycle, ewcMoveE);
858 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
859 inputrec, nullptr, nullptr, nullptr, 1, &terminate,
860 nullptr, FALSE,
861 CGLO_ENERGY |
862 CGLO_PRESSURE |
863 CGLO_CONSTRAINT);
865 wallcycle_stop(wcycle, ewcMoveE);
868 /* Calculate long range corrections to pressure and energy */
869 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
870 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
871 enerd->term[F_DISPCORR] = enercorr;
872 enerd->term[F_EPOT] += enercorr;
873 enerd->term[F_PRES] += prescorr;
874 enerd->term[F_DVDL] += dvdlcorr;
876 ems->epot = enerd->term[F_EPOT];
878 if (constr)
880 /* Project out the constraint components of the force */
881 wallcycle_start(wcycle, ewcCONSTR);
882 dvdl_constr = 0;
883 rvec *f_rvec = as_rvec_array(ems->f.data());
884 constrain(nullptr, FALSE, FALSE, constr, &top->idef,
885 inputrec, cr, ms, count, 0, 1.0, mdAtoms->mdatoms(),
886 as_rvec_array(ems->s.x.data()), f_rvec, f_rvec,
887 fr->bMolPBC, ems->s.box,
888 ems->s.lambda[efptBONDED], &dvdl_constr,
889 nullptr, &shake_vir, nrnb, gmx::econqForceDispl);
890 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
891 m_add(force_vir, shake_vir, vir);
892 wallcycle_stop(wcycle, ewcCONSTR);
894 else
896 copy_mat(force_vir, vir);
899 clear_mat(ekin);
900 enerd->term[F_PRES] =
901 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
903 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
905 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
907 get_state_f_norm_max(cr, &(inputrec->opts), mdAtoms->mdatoms(), ems);
911 } // namespace
913 //! Parallel utility summing energies and forces
914 static double reorder_partsum(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
915 gmx_mtop_t *top_global,
916 em_state_t *s_min, em_state_t *s_b)
918 t_block *cgs_gl;
919 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
920 double partsum;
921 unsigned char *grpnrFREEZE;
923 if (debug)
925 fprintf(debug, "Doing reorder_partsum\n");
928 const rvec *fm = as_rvec_array(s_min->f.data());
929 const rvec *fb = as_rvec_array(s_b->f.data());
931 cgs_gl = dd_charge_groups_global(cr->dd);
932 index = cgs_gl->index;
934 /* Collect fm in a global vector fmg.
935 * This conflicts with the spirit of domain decomposition,
936 * but to fully optimize this a much more complicated algorithm is required.
938 rvec *fmg;
939 snew(fmg, top_global->natoms);
941 ncg = s_min->s.cg_gl.size();
942 cg_gl = s_min->s.cg_gl.data();
943 i = 0;
944 for (c = 0; c < ncg; c++)
946 cg = cg_gl[c];
947 a0 = index[cg];
948 a1 = index[cg+1];
949 for (a = a0; a < a1; a++)
951 copy_rvec(fm[i], fmg[a]);
952 i++;
955 gmx_sum(top_global->natoms*3, fmg[0], cr);
957 /* Now we will determine the part of the sum for the cgs in state s_b */
958 ncg = s_b->s.cg_gl.size();
959 cg_gl = s_b->s.cg_gl.data();
960 partsum = 0;
961 i = 0;
962 gf = 0;
963 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
964 for (c = 0; c < ncg; c++)
966 cg = cg_gl[c];
967 a0 = index[cg];
968 a1 = index[cg+1];
969 for (a = a0; a < a1; a++)
971 if (mdatoms->cFREEZE && grpnrFREEZE)
973 gf = grpnrFREEZE[i];
975 for (m = 0; m < DIM; m++)
977 if (!opts->nFreeze[gf][m])
979 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
982 i++;
986 sfree(fmg);
988 return partsum;
991 //! Print some stuff, like beta, whatever that means.
992 static real pr_beta(const t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
993 gmx_mtop_t *top_global,
994 em_state_t *s_min, em_state_t *s_b)
996 double sum;
998 /* This is just the classical Polak-Ribiere calculation of beta;
999 * it looks a bit complicated since we take freeze groups into account,
1000 * and might have to sum it in parallel runs.
1003 if (!DOMAINDECOMP(cr) ||
1004 (s_min->s.ddp_count == cr->dd->ddp_count &&
1005 s_b->s.ddp_count == cr->dd->ddp_count))
1007 const rvec *fm = as_rvec_array(s_min->f.data());
1008 const rvec *fb = as_rvec_array(s_b->f.data());
1009 sum = 0;
1010 int gf = 0;
1011 /* This part of code can be incorrect with DD,
1012 * since the atom ordering in s_b and s_min might differ.
1014 for (int i = 0; i < mdatoms->homenr; i++)
1016 if (mdatoms->cFREEZE)
1018 gf = mdatoms->cFREEZE[i];
1020 for (int m = 0; m < DIM; m++)
1022 if (!opts->nFreeze[gf][m])
1024 sum += (fb[i][m] - fm[i][m])*fb[i][m];
1029 else
1031 /* We need to reorder cgs while summing */
1032 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
1034 if (PAR(cr))
1036 gmx_sumd(1, &sum, cr);
1039 return sum/gmx::square(s_min->fnorm);
1042 namespace gmx
1045 void
1046 Integrator::do_cg()
1048 const char *CG = "Polak-Ribiere Conjugate Gradients";
1050 gmx_localtop_t *top;
1051 gmx_enerdata_t *enerd;
1052 gmx_global_stat_t gstat;
1053 t_graph *graph;
1054 double tmp, minstep;
1055 real stepsize;
1056 real a, b, c, beta = 0.0;
1057 real epot_repl = 0;
1058 real pnorm;
1059 t_mdebin *mdebin;
1060 gmx_bool converged, foundlower;
1061 rvec mu_tot;
1062 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1063 tensor vir, pres;
1064 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1065 gmx_mdoutf_t outf;
1066 int m, step, nminstep;
1067 auto mdatoms = mdAtoms->mdatoms();
1069 step = 0;
1071 // Ensure the extra per-atom state array gets allocated
1072 state_global->flags |= (1<<estCGP);
1074 /* Create 4 states on the stack and extract pointers that we will swap */
1075 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1076 em_state_t *s_min = &s0;
1077 em_state_t *s_a = &s1;
1078 em_state_t *s_b = &s2;
1079 em_state_t *s_c = &s3;
1081 /* Init em and store the local state in s_min */
1082 init_em(fplog, CG, cr, ms, outputProvider, inputrec, mdrunOptions,
1083 state_global, top_global, s_min, &top,
1084 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1085 vsite, constr, nullptr,
1086 nfile, fnm, &outf, &mdebin, wcycle);
1088 /* Print to log file */
1089 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1091 /* Max number of steps */
1092 number_steps = inputrec->nsteps;
1094 if (MASTER(cr))
1096 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1098 if (fplog)
1100 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1103 EnergyEvaluator energyEvaluator {
1104 fplog, cr, ms,
1105 top_global, top,
1106 inputrec, nrnb, wcycle, gstat,
1107 vsite, constr, fcd, graph,
1108 mdAtoms, fr, enerd
1110 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1111 /* do_force always puts the charge groups in the box and shifts again
1112 * We do not unshift, so molecules are always whole in congrad.c
1114 energyEvaluator.run(s_min, mu_tot, vir, pres, -1, TRUE);
1116 if (MASTER(cr))
1118 /* Copy stuff to the energy bin for easy printing etc. */
1119 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1120 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1121 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1123 print_ebin_header(fplog, step, step);
1124 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1125 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1128 /* Estimate/guess the initial stepsize */
1129 stepsize = inputrec->em_stepsize/s_min->fnorm;
1131 if (MASTER(cr))
1133 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1134 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1135 s_min->fmax, s_min->a_fmax+1);
1136 fprintf(stderr, " F-Norm = %12.5e\n",
1137 s_min->fnorm/sqrtNumAtoms);
1138 fprintf(stderr, "\n");
1139 /* and copy to the log file too... */
1140 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1141 s_min->fmax, s_min->a_fmax+1);
1142 fprintf(fplog, " F-Norm = %12.5e\n",
1143 s_min->fnorm/sqrtNumAtoms);
1144 fprintf(fplog, "\n");
1146 /* Start the loop over CG steps.
1147 * Each successful step is counted, and we continue until
1148 * we either converge or reach the max number of steps.
1150 converged = FALSE;
1151 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1154 /* start taking steps in a new direction
1155 * First time we enter the routine, beta=0, and the direction is
1156 * simply the negative gradient.
1159 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1160 rvec *pm = as_rvec_array(s_min->s.cg_p.data());
1161 const rvec *sfm = as_rvec_array(s_min->f.data());
1162 double gpa = 0;
1163 int gf = 0;
1164 for (int i = 0; i < mdatoms->homenr; i++)
1166 if (mdatoms->cFREEZE)
1168 gf = mdatoms->cFREEZE[i];
1170 for (m = 0; m < DIM; m++)
1172 if (!inputrec->opts.nFreeze[gf][m])
1174 pm[i][m] = sfm[i][m] + beta*pm[i][m];
1175 gpa -= pm[i][m]*sfm[i][m];
1176 /* f is negative gradient, thus the sign */
1178 else
1180 pm[i][m] = 0;
1185 /* Sum the gradient along the line across CPUs */
1186 if (PAR(cr))
1188 gmx_sumd(1, &gpa, cr);
1191 /* Calculate the norm of the search vector */
1192 get_f_norm_max(cr, &(inputrec->opts), mdatoms, pm, &pnorm, nullptr, nullptr);
1194 /* Just in case stepsize reaches zero due to numerical precision... */
1195 if (stepsize <= 0)
1197 stepsize = inputrec->em_stepsize/pnorm;
1201 * Double check the value of the derivative in the search direction.
1202 * If it is positive it must be due to the old information in the
1203 * CG formula, so just remove that and start over with beta=0.
1204 * This corresponds to a steepest descent step.
1206 if (gpa > 0)
1208 beta = 0;
1209 step--; /* Don't count this step since we are restarting */
1210 continue; /* Go back to the beginning of the big for-loop */
1213 /* Calculate minimum allowed stepsize, before the average (norm)
1214 * relative change in coordinate is smaller than precision
1216 minstep = 0;
1217 for (int i = 0; i < mdatoms->homenr; i++)
1219 for (m = 0; m < DIM; m++)
1221 tmp = fabs(s_min->s.x[i][m]);
1222 if (tmp < 1.0)
1224 tmp = 1.0;
1226 tmp = pm[i][m]/tmp;
1227 minstep += tmp*tmp;
1230 /* Add up from all CPUs */
1231 if (PAR(cr))
1233 gmx_sumd(1, &minstep, cr);
1236 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1238 if (stepsize < minstep)
1240 converged = TRUE;
1241 break;
1244 /* Write coordinates if necessary */
1245 do_x = do_per_step(step, inputrec->nstxout);
1246 do_f = do_per_step(step, inputrec->nstfout);
1248 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
1249 top_global, inputrec, step,
1250 s_min, state_global, observablesHistory);
1252 /* Take a step downhill.
1253 * In theory, we should minimize the function along this direction.
1254 * That is quite possible, but it turns out to take 5-10 function evaluations
1255 * for each line. However, we dont really need to find the exact minimum -
1256 * it is much better to start a new CG step in a modified direction as soon
1257 * as we are close to it. This will save a lot of energy evaluations.
1259 * In practice, we just try to take a single step.
1260 * If it worked (i.e. lowered the energy), we increase the stepsize but
1261 * the continue straight to the next CG step without trying to find any minimum.
1262 * If it didn't work (higher energy), there must be a minimum somewhere between
1263 * the old position and the new one.
1265 * Due to the finite numerical accuracy, it turns out that it is a good idea
1266 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1267 * This leads to lower final energies in the tests I've done. / Erik
1269 s_a->epot = s_min->epot;
1270 a = 0.0;
1271 c = a + stepsize; /* reference position along line is zero */
1273 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1275 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1276 s_min, top, mdAtoms, fr, vsite, constr,
1277 nrnb, wcycle);
1280 /* Take a trial step (new coords in s_c) */
1281 do_em_step(cr, ms, inputrec, mdatoms, fr->bMolPBC, s_min, c, &s_min->s.cg_p, s_c,
1282 constr, top, nrnb, wcycle, -1);
1284 neval++;
1285 /* Calculate energy for the trial step */
1286 energyEvaluator.run(s_c, mu_tot, vir, pres, -1, FALSE);
1288 /* Calc derivative along line */
1289 const rvec *pc = as_rvec_array(s_c->s.cg_p.data());
1290 const rvec *sfc = as_rvec_array(s_c->f.data());
1291 double gpc = 0;
1292 for (int i = 0; i < mdatoms->homenr; i++)
1294 for (m = 0; m < DIM; m++)
1296 gpc -= pc[i][m]*sfc[i][m]; /* f is negative gradient, thus the sign */
1299 /* Sum the gradient along the line across CPUs */
1300 if (PAR(cr))
1302 gmx_sumd(1, &gpc, cr);
1305 /* This is the max amount of increase in energy we tolerate */
1306 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1308 /* Accept the step if the energy is lower, or if it is not significantly higher
1309 * and the line derivative is still negative.
1311 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1313 foundlower = TRUE;
1314 /* Great, we found a better energy. Increase step for next iteration
1315 * if we are still going down, decrease it otherwise
1317 if (gpc < 0)
1319 stepsize *= 1.618034; /* The golden section */
1321 else
1323 stepsize *= 0.618034; /* 1/golden section */
1326 else
1328 /* New energy is the same or higher. We will have to do some work
1329 * to find a smaller value in the interval. Take smaller step next time!
1331 foundlower = FALSE;
1332 stepsize *= 0.618034;
1338 /* OK, if we didn't find a lower value we will have to locate one now - there must
1339 * be one in the interval [a=0,c].
1340 * The same thing is valid here, though: Don't spend dozens of iterations to find
1341 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1342 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1344 * I also have a safeguard for potentially really pathological functions so we never
1345 * take more than 20 steps before we give up ...
1347 * If we already found a lower value we just skip this step and continue to the update.
1349 double gpb;
1350 if (!foundlower)
1352 nminstep = 0;
1356 /* Select a new trial point.
1357 * If the derivatives at points a & c have different sign we interpolate to zero,
1358 * otherwise just do a bisection.
1360 if (gpa < 0 && gpc > 0)
1362 b = a + gpa*(a-c)/(gpc-gpa);
1364 else
1366 b = 0.5*(a+c);
1369 /* safeguard if interpolation close to machine accuracy causes errors:
1370 * never go outside the interval
1372 if (b <= a || b >= c)
1374 b = 0.5*(a+c);
1377 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1379 /* Reload the old state */
1380 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1381 s_min, top, mdAtoms, fr, vsite, constr,
1382 nrnb, wcycle);
1385 /* Take a trial step to this new point - new coords in s_b */
1386 do_em_step(cr, ms, inputrec, mdatoms, fr->bMolPBC, s_min, b, &s_min->s.cg_p, s_b,
1387 constr, top, nrnb, wcycle, -1);
1389 neval++;
1390 /* Calculate energy for the trial step */
1391 energyEvaluator.run(s_b, mu_tot, vir, pres, -1, FALSE);
1393 /* p does not change within a step, but since the domain decomposition
1394 * might change, we have to use cg_p of s_b here.
1396 const rvec *pb = as_rvec_array(s_b->s.cg_p.data());
1397 const rvec *sfb = as_rvec_array(s_b->f.data());
1398 gpb = 0;
1399 for (int i = 0; i < mdatoms->homenr; i++)
1401 for (m = 0; m < DIM; m++)
1403 gpb -= pb[i][m]*sfb[i][m]; /* f is negative gradient, thus the sign */
1406 /* Sum the gradient along the line across CPUs */
1407 if (PAR(cr))
1409 gmx_sumd(1, &gpb, cr);
1412 if (debug)
1414 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1415 s_a->epot, s_b->epot, s_c->epot, gpb);
1418 epot_repl = s_b->epot;
1420 /* Keep one of the intervals based on the value of the derivative at the new point */
1421 if (gpb > 0)
1423 /* Replace c endpoint with b */
1424 swap_em_state(&s_b, &s_c);
1425 c = b;
1426 gpc = gpb;
1428 else
1430 /* Replace a endpoint with b */
1431 swap_em_state(&s_b, &s_a);
1432 a = b;
1433 gpa = gpb;
1437 * Stop search as soon as we find a value smaller than the endpoints.
1438 * Never run more than 20 steps, no matter what.
1440 nminstep++;
1442 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1443 (nminstep < 20));
1445 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1446 nminstep >= 20)
1448 /* OK. We couldn't find a significantly lower energy.
1449 * If beta==0 this was steepest descent, and then we give up.
1450 * If not, set beta=0 and restart with steepest descent before quitting.
1452 if (beta == 0.0)
1454 /* Converged */
1455 converged = TRUE;
1456 break;
1458 else
1460 /* Reset memory before giving up */
1461 beta = 0.0;
1462 continue;
1466 /* Select min energy state of A & C, put the best in B.
1468 if (s_c->epot < s_a->epot)
1470 if (debug)
1472 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1473 s_c->epot, s_a->epot);
1475 swap_em_state(&s_b, &s_c);
1476 gpb = gpc;
1478 else
1480 if (debug)
1482 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1483 s_a->epot, s_c->epot);
1485 swap_em_state(&s_b, &s_a);
1486 gpb = gpa;
1490 else
1492 if (debug)
1494 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1495 s_c->epot);
1497 swap_em_state(&s_b, &s_c);
1498 gpb = gpc;
1501 /* new search direction */
1502 /* beta = 0 means forget all memory and restart with steepest descents. */
1503 if (nstcg && ((step % nstcg) == 0))
1505 beta = 0.0;
1507 else
1509 /* s_min->fnorm cannot be zero, because then we would have converged
1510 * and broken out.
1513 /* Polak-Ribiere update.
1514 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1516 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1518 /* Limit beta to prevent oscillations */
1519 if (fabs(beta) > 5.0)
1521 beta = 0.0;
1525 /* update positions */
1526 swap_em_state(&s_min, &s_b);
1527 gpa = gpb;
1529 /* Print it if necessary */
1530 if (MASTER(cr))
1532 if (mdrunOptions.verbose)
1534 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1535 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1536 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1537 s_min->fmax, s_min->a_fmax+1);
1538 fflush(stderr);
1540 /* Store the new (lower) energies */
1541 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1542 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1543 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1545 do_log = do_per_step(step, inputrec->nstlog);
1546 do_ene = do_per_step(step, inputrec->nstenergy);
1548 /* Prepare IMD energy record, if bIMD is TRUE. */
1549 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1551 if (do_log)
1553 print_ebin_header(fplog, step, step);
1555 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1556 do_log ? fplog : nullptr, step, step, eprNORMAL,
1557 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1560 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1561 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
1563 IMD_send_positions(inputrec->imd);
1566 /* Stop when the maximum force lies below tolerance.
1567 * If we have reached machine precision, converged is already set to true.
1569 converged = converged || (s_min->fmax < inputrec->em_tol);
1571 } /* End of the loop */
1573 /* IMD cleanup, if bIMD is TRUE. */
1574 IMD_finalize(inputrec->bIMD, inputrec->imd);
1576 if (converged)
1578 step--; /* we never took that last step in this case */
1581 if (s_min->fmax > inputrec->em_tol)
1583 if (MASTER(cr))
1585 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1586 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1588 converged = FALSE;
1591 if (MASTER(cr))
1593 /* If we printed energy and/or logfile last step (which was the last step)
1594 * we don't have to do it again, but otherwise print the final values.
1596 if (!do_log)
1598 /* Write final value to log since we didn't do anything the last step */
1599 print_ebin_header(fplog, step, step);
1601 if (!do_ene || !do_log)
1603 /* Write final energy file entries */
1604 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1605 !do_log ? fplog : nullptr, step, step, eprNORMAL,
1606 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1610 /* Print some stuff... */
1611 if (MASTER(cr))
1613 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1616 /* IMPORTANT!
1617 * For accurate normal mode calculation it is imperative that we
1618 * store the last conformation into the full precision binary trajectory.
1620 * However, we should only do it if we did NOT already write this step
1621 * above (which we did if do_x or do_f was true).
1623 do_x = !do_per_step(step, inputrec->nstxout);
1624 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1626 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1627 top_global, inputrec, step,
1628 s_min, state_global, observablesHistory);
1631 if (MASTER(cr))
1633 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1634 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1635 s_min, sqrtNumAtoms);
1636 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1637 s_min, sqrtNumAtoms);
1639 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1642 finish_em(cr, outf, walltime_accounting, wcycle);
1644 /* To print the actual number of steps we needed somewhere */
1645 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1649 void
1650 Integrator::do_lbfgs()
1652 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1653 em_state_t ems;
1654 gmx_localtop_t *top;
1655 gmx_enerdata_t *enerd;
1656 gmx_global_stat_t gstat;
1657 t_graph *graph;
1658 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1659 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1660 real *rho, *alpha, *p, *s, **dx, **dg;
1661 real a, b, c, maxdelta, delta;
1662 real diag, Epot0;
1663 real dgdx, dgdg, sq, yr, beta;
1664 t_mdebin *mdebin;
1665 gmx_bool converged;
1666 rvec mu_tot;
1667 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1668 tensor vir, pres;
1669 int start, end, number_steps;
1670 gmx_mdoutf_t outf;
1671 int i, k, m, n, gf, step;
1672 int mdof_flags;
1673 auto mdatoms = mdAtoms->mdatoms();
1675 if (PAR(cr))
1677 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1680 if (nullptr != constr)
1682 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).");
1685 n = 3*state_global->natoms;
1686 nmaxcorr = inputrec->nbfgscorr;
1688 snew(frozen, n);
1690 snew(p, n);
1691 snew(rho, nmaxcorr);
1692 snew(alpha, nmaxcorr);
1694 snew(dx, nmaxcorr);
1695 for (i = 0; i < nmaxcorr; i++)
1697 snew(dx[i], n);
1700 snew(dg, nmaxcorr);
1701 for (i = 0; i < nmaxcorr; i++)
1703 snew(dg[i], n);
1706 step = 0;
1707 neval = 0;
1709 /* Init em */
1710 init_em(fplog, LBFGS, cr, ms, outputProvider, inputrec, mdrunOptions,
1711 state_global, top_global, &ems, &top,
1712 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
1713 vsite, constr, nullptr,
1714 nfile, fnm, &outf, &mdebin, wcycle);
1716 start = 0;
1717 end = mdatoms->homenr;
1719 /* We need 4 working states */
1720 em_state_t s0 {}, s1 {}, s2 {}, s3 {};
1721 em_state_t *sa = &s0;
1722 em_state_t *sb = &s1;
1723 em_state_t *sc = &s2;
1724 em_state_t *last = &s3;
1725 /* Initialize by copying the state from ems (we could skip x and f here) */
1726 *sa = ems;
1727 *sb = ems;
1728 *sc = ems;
1730 /* Print to log file */
1731 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1733 do_log = do_ene = do_x = do_f = TRUE;
1735 /* Max number of steps */
1736 number_steps = inputrec->nsteps;
1738 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1739 gf = 0;
1740 for (i = start; i < end; i++)
1742 if (mdatoms->cFREEZE)
1744 gf = mdatoms->cFREEZE[i];
1746 for (m = 0; m < DIM; m++)
1748 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1751 if (MASTER(cr))
1753 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1755 if (fplog)
1757 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1760 if (vsite)
1762 construct_vsites(vsite, as_rvec_array(state_global->x.data()), 1, nullptr,
1763 top->idef.iparams, top->idef.il,
1764 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1767 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1768 /* do_force always puts the charge groups in the box and shifts again
1769 * We do not unshift, so molecules are always whole
1771 neval++;
1772 EnergyEvaluator energyEvaluator {
1773 fplog, cr, ms,
1774 top_global, top,
1775 inputrec, nrnb, wcycle, gstat,
1776 vsite, constr, fcd, graph,
1777 mdAtoms, fr, enerd
1779 energyEvaluator.run(&ems, mu_tot, vir, pres, -1, TRUE);
1781 if (MASTER(cr))
1783 /* Copy stuff to the energy bin for easy printing etc. */
1784 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1785 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
1786 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
1788 print_ebin_header(fplog, step, step);
1789 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1790 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
1793 /* Set the initial step.
1794 * since it will be multiplied by the non-normalized search direction
1795 * vector (force vector the first time), we scale it by the
1796 * norm of the force.
1799 if (MASTER(cr))
1801 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1802 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1803 fprintf(stderr, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1804 fprintf(stderr, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1805 fprintf(stderr, "\n");
1806 /* and copy to the log file too... */
1807 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1808 fprintf(fplog, " F-max = %12.5e on atom %d\n", ems.fmax, ems.a_fmax + 1);
1809 fprintf(fplog, " F-Norm = %12.5e\n", ems.fnorm/sqrtNumAtoms);
1810 fprintf(fplog, "\n");
1813 // Point is an index to the memory of search directions, where 0 is the first one.
1814 point = 0;
1816 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1817 real *fInit = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1818 for (i = 0; i < n; i++)
1820 if (!frozen[i])
1822 dx[point][i] = fInit[i]; /* Initial search direction */
1824 else
1826 dx[point][i] = 0;
1830 // Stepsize will be modified during the search, and actually it is not critical
1831 // (the main efficiency in the algorithm comes from changing directions), but
1832 // we still need an initial value, so estimate it as the inverse of the norm
1833 // so we take small steps where the potential fluctuates a lot.
1834 stepsize = 1.0/ems.fnorm;
1836 /* Start the loop over BFGS steps.
1837 * Each successful step is counted, and we continue until
1838 * we either converge or reach the max number of steps.
1841 ncorr = 0;
1843 /* Set the gradient from the force */
1844 converged = FALSE;
1845 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1848 /* Write coordinates if necessary */
1849 do_x = do_per_step(step, inputrec->nstxout);
1850 do_f = do_per_step(step, inputrec->nstfout);
1852 mdof_flags = 0;
1853 if (do_x)
1855 mdof_flags |= MDOF_X;
1858 if (do_f)
1860 mdof_flags |= MDOF_F;
1863 if (inputrec->bIMD)
1865 mdof_flags |= MDOF_IMD;
1868 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1869 top_global, step, (real)step, &ems.s, state_global, observablesHistory, ems.f);
1871 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1873 /* make s a pointer to current search direction - point=0 first time we get here */
1874 s = dx[point];
1876 real *xx = static_cast<real *>(as_rvec_array(ems.s.x.data())[0]);
1877 real *ff = static_cast<real *>(as_rvec_array(ems.f.data())[0]);
1879 // calculate line gradient in position A
1880 for (gpa = 0, i = 0; i < n; i++)
1882 gpa -= s[i]*ff[i];
1885 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1886 * relative change in coordinate is smaller than precision
1888 for (minstep = 0, i = 0; i < n; i++)
1890 tmp = fabs(xx[i]);
1891 if (tmp < 1.0)
1893 tmp = 1.0;
1895 tmp = s[i]/tmp;
1896 minstep += tmp*tmp;
1898 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1900 if (stepsize < minstep)
1902 converged = TRUE;
1903 break;
1906 // Before taking any steps along the line, store the old position
1907 *last = ems;
1908 real *lastx = static_cast<real *>(as_rvec_array(last->s.x.data())[0]);
1909 real *lastf = static_cast<real *>(as_rvec_array(last->f.data())[0]);
1910 Epot0 = ems.epot;
1912 *sa = ems;
1914 /* Take a step downhill.
1915 * In theory, we should find the actual minimum of the function in this
1916 * direction, somewhere along the line.
1917 * That is quite possible, but it turns out to take 5-10 function evaluations
1918 * for each line. However, we dont really need to find the exact minimum -
1919 * it is much better to start a new BFGS step in a modified direction as soon
1920 * as we are close to it. This will save a lot of energy evaluations.
1922 * In practice, we just try to take a single step.
1923 * If it worked (i.e. lowered the energy), we increase the stepsize but
1924 * continue straight to the next BFGS step without trying to find any minimum,
1925 * i.e. we change the search direction too. If the line was smooth, it is
1926 * likely we are in a smooth region, and then it makes sense to take longer
1927 * steps in the modified search direction too.
1929 * If it didn't work (higher energy), there must be a minimum somewhere between
1930 * the old position and the new one. Then we need to start by finding a lower
1931 * value before we change search direction. Since the energy was apparently
1932 * quite rough, we need to decrease the step size.
1934 * Due to the finite numerical accuracy, it turns out that it is a good idea
1935 * to accept a SMALL increase in energy, if the derivative is still downhill.
1936 * This leads to lower final energies in the tests I've done. / Erik
1939 // State "A" is the first position along the line.
1940 // reference position along line is initially zero
1941 a = 0.0;
1943 // Check stepsize first. We do not allow displacements
1944 // larger than emstep.
1948 // Pick a new position C by adding stepsize to A.
1949 c = a + stepsize;
1951 // Calculate what the largest change in any individual coordinate
1952 // would be (translation along line * gradient along line)
1953 maxdelta = 0;
1954 for (i = 0; i < n; i++)
1956 delta = c*s[i];
1957 if (delta > maxdelta)
1959 maxdelta = delta;
1962 // If any displacement is larger than the stepsize limit, reduce the step
1963 if (maxdelta > inputrec->em_stepsize)
1965 stepsize *= 0.1;
1968 while (maxdelta > inputrec->em_stepsize);
1970 // Take a trial step and move the coordinate array xc[] to position C
1971 real *xc = static_cast<real *>(as_rvec_array(sc->s.x.data())[0]);
1972 for (i = 0; i < n; i++)
1974 xc[i] = lastx[i] + c*s[i];
1977 neval++;
1978 // Calculate energy for the trial step in position C
1979 energyEvaluator.run(sc, mu_tot, vir, pres, step, FALSE);
1981 // Calc line gradient in position C
1982 real *fc = static_cast<real *>(as_rvec_array(sc->f.data())[0]);
1983 for (gpc = 0, i = 0; i < n; i++)
1985 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1987 /* Sum the gradient along the line across CPUs */
1988 if (PAR(cr))
1990 gmx_sumd(1, &gpc, cr);
1993 // This is the max amount of increase in energy we tolerate.
1994 // By allowing VERY small changes (close to numerical precision) we
1995 // frequently find even better (lower) final energies.
1996 tmp = sqrt(GMX_REAL_EPS)*fabs(sa->epot);
1998 // Accept the step if the energy is lower in the new position C (compared to A),
1999 // or if it is not significantly higher and the line derivative is still negative.
2000 if (sc->epot < sa->epot || (gpc < 0 && sc->epot < (sa->epot + tmp)))
2002 // Great, we found a better energy. We no longer try to alter the
2003 // stepsize, but simply accept this new better position. The we select a new
2004 // search direction instead, which will be much more efficient than continuing
2005 // to take smaller steps along a line. Set fnorm based on the new C position,
2006 // which will be used to update the stepsize to 1/fnorm further down.
2007 foundlower = TRUE;
2009 else
2011 // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2012 // or higher than in point A. In this case it is pointless to move to point C,
2013 // so we will have to do more iterations along the same line to find a smaller
2014 // value in the interval [A=0.0,C].
2015 // Here, A is still 0.0, but that will change when we do a search in the interval
2016 // [0.0,C] below. That search we will do by interpolation or bisection rather
2017 // than with the stepsize, so no need to modify it. For the next search direction
2018 // it will be reset to 1/fnorm anyway.
2019 foundlower = FALSE;
2022 if (!foundlower)
2024 // OK, if we didn't find a lower value we will have to locate one now - there must
2025 // be one in the interval [a,c].
2026 // The same thing is valid here, though: Don't spend dozens of iterations to find
2027 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2028 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2029 // I also have a safeguard for potentially really pathological functions so we never
2030 // take more than 20 steps before we give up.
2031 // If we already found a lower value we just skip this step and continue to the update.
2032 real fnorm = 0;
2033 nminstep = 0;
2036 // Select a new trial point B in the interval [A,C].
2037 // If the derivatives at points a & c have different sign we interpolate to zero,
2038 // otherwise just do a bisection since there might be multiple minima/maxima
2039 // inside the interval.
2040 if (gpa < 0 && gpc > 0)
2042 b = a + gpa*(a-c)/(gpc-gpa);
2044 else
2046 b = 0.5*(a+c);
2049 /* safeguard if interpolation close to machine accuracy causes errors:
2050 * never go outside the interval
2052 if (b <= a || b >= c)
2054 b = 0.5*(a+c);
2057 // Take a trial step to point B
2058 real *xb = static_cast<real *>(as_rvec_array(sb->s.x.data())[0]);
2059 for (i = 0; i < n; i++)
2061 xb[i] = lastx[i] + b*s[i];
2064 neval++;
2065 // Calculate energy for the trial step in point B
2066 energyEvaluator.run(sb, mu_tot, vir, pres, step, FALSE);
2067 fnorm = sb->fnorm;
2069 // Calculate gradient in point B
2070 real *fb = static_cast<real *>(as_rvec_array(sb->f.data())[0]);
2071 for (gpb = 0, i = 0; i < n; i++)
2073 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2076 /* Sum the gradient along the line across CPUs */
2077 if (PAR(cr))
2079 gmx_sumd(1, &gpb, cr);
2082 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2083 // at the new point B, and rename the endpoints of this new interval A and C.
2084 if (gpb > 0)
2086 /* Replace c endpoint with b */
2087 c = b;
2088 /* swap states b and c */
2089 swap_em_state(&sb, &sc);
2091 else
2093 /* Replace a endpoint with b */
2094 a = b;
2095 /* swap states a and b */
2096 swap_em_state(&sa, &sb);
2100 * Stop search as soon as we find a value smaller than the endpoints,
2101 * or if the tolerance is below machine precision.
2102 * Never run more than 20 steps, no matter what.
2104 nminstep++;
2106 while ((sb->epot > sa->epot || sb->epot > sc->epot) && (nminstep < 20));
2108 if (fabs(sb->epot - Epot0) < GMX_REAL_EPS || nminstep >= 20)
2110 /* OK. We couldn't find a significantly lower energy.
2111 * If ncorr==0 this was steepest descent, and then we give up.
2112 * If not, reset memory to restart as steepest descent before quitting.
2114 if (ncorr == 0)
2116 /* Converged */
2117 converged = TRUE;
2118 break;
2120 else
2122 /* Reset memory */
2123 ncorr = 0;
2124 /* Search in gradient direction */
2125 for (i = 0; i < n; i++)
2127 dx[point][i] = ff[i];
2129 /* Reset stepsize */
2130 stepsize = 1.0/fnorm;
2131 continue;
2135 /* Select min energy state of A & C, put the best in xx/ff/Epot
2137 if (sc->epot < sa->epot)
2139 /* Use state C */
2140 ems = *sc;
2141 step_taken = c;
2143 else
2145 /* Use state A */
2146 ems = *sa;
2147 step_taken = a;
2151 else
2153 /* found lower */
2154 /* Use state C */
2155 ems = *sc;
2156 step_taken = c;
2159 /* Update the memory information, and calculate a new
2160 * approximation of the inverse hessian
2163 /* Have new data in Epot, xx, ff */
2164 if (ncorr < nmaxcorr)
2166 ncorr++;
2169 for (i = 0; i < n; i++)
2171 dg[point][i] = lastf[i]-ff[i];
2172 dx[point][i] *= step_taken;
2175 dgdg = 0;
2176 dgdx = 0;
2177 for (i = 0; i < n; i++)
2179 dgdg += dg[point][i]*dg[point][i];
2180 dgdx += dg[point][i]*dx[point][i];
2183 diag = dgdx/dgdg;
2185 rho[point] = 1.0/dgdx;
2186 point++;
2188 if (point >= nmaxcorr)
2190 point = 0;
2193 /* Update */
2194 for (i = 0; i < n; i++)
2196 p[i] = ff[i];
2199 cp = point;
2201 /* Recursive update. First go back over the memory points */
2202 for (k = 0; k < ncorr; k++)
2204 cp--;
2205 if (cp < 0)
2207 cp = ncorr-1;
2210 sq = 0;
2211 for (i = 0; i < n; i++)
2213 sq += dx[cp][i]*p[i];
2216 alpha[cp] = rho[cp]*sq;
2218 for (i = 0; i < n; i++)
2220 p[i] -= alpha[cp]*dg[cp][i];
2224 for (i = 0; i < n; i++)
2226 p[i] *= diag;
2229 /* And then go forward again */
2230 for (k = 0; k < ncorr; k++)
2232 yr = 0;
2233 for (i = 0; i < n; i++)
2235 yr += p[i]*dg[cp][i];
2238 beta = rho[cp]*yr;
2239 beta = alpha[cp]-beta;
2241 for (i = 0; i < n; i++)
2243 p[i] += beta*dx[cp][i];
2246 cp++;
2247 if (cp >= ncorr)
2249 cp = 0;
2253 for (i = 0; i < n; i++)
2255 if (!frozen[i])
2257 dx[point][i] = p[i];
2259 else
2261 dx[point][i] = 0;
2265 /* Print it if necessary */
2266 if (MASTER(cr))
2268 if (mdrunOptions.verbose)
2270 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2271 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2272 step, ems.epot, ems.fnorm/sqrtNumAtoms, ems.fmax, ems.a_fmax + 1);
2273 fflush(stderr);
2275 /* Store the new (lower) energies */
2276 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2277 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
2278 nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2279 do_log = do_per_step(step, inputrec->nstlog);
2280 do_ene = do_per_step(step, inputrec->nstenergy);
2281 if (do_log)
2283 print_ebin_header(fplog, step, step);
2285 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2286 do_log ? fplog : nullptr, step, step, eprNORMAL,
2287 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2290 /* Send x and E to IMD client, if bIMD is TRUE. */
2291 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, as_rvec_array(state_global->x.data()), inputrec, 0, wcycle) && MASTER(cr))
2293 IMD_send_positions(inputrec->imd);
2296 // Reset stepsize in we are doing more iterations
2297 stepsize = 1.0/ems.fnorm;
2299 /* Stop when the maximum force lies below tolerance.
2300 * If we have reached machine precision, converged is already set to true.
2302 converged = converged || (ems.fmax < inputrec->em_tol);
2304 } /* End of the loop */
2306 /* IMD cleanup, if bIMD is TRUE. */
2307 IMD_finalize(inputrec->bIMD, inputrec->imd);
2309 if (converged)
2311 step--; /* we never took that last step in this case */
2314 if (ems.fmax > inputrec->em_tol)
2316 if (MASTER(cr))
2318 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2319 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2321 converged = FALSE;
2324 /* If we printed energy and/or logfile last step (which was the last step)
2325 * we don't have to do it again, but otherwise print the final values.
2327 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2329 print_ebin_header(fplog, step, step);
2331 if (!do_ene || !do_log) /* Write final energy file entries */
2333 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2334 !do_log ? fplog : nullptr, step, step, eprNORMAL,
2335 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2338 /* Print some stuff... */
2339 if (MASTER(cr))
2341 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2344 /* IMPORTANT!
2345 * For accurate normal mode calculation it is imperative that we
2346 * store the last conformation into the full precision binary trajectory.
2348 * However, we should only do it if we did NOT already write this step
2349 * above (which we did if do_x or do_f was true).
2351 do_x = !do_per_step(step, inputrec->nstxout);
2352 do_f = !do_per_step(step, inputrec->nstfout);
2353 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2354 top_global, inputrec, step,
2355 &ems, state_global, observablesHistory);
2357 if (MASTER(cr))
2359 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2360 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2361 number_steps, &ems, sqrtNumAtoms);
2362 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2363 number_steps, &ems, sqrtNumAtoms);
2365 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2368 finish_em(cr, outf, walltime_accounting, wcycle);
2370 /* To print the actual number of steps we needed somewhere */
2371 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2374 void
2375 Integrator::do_steep()
2377 const char *SD = "Steepest Descents";
2378 gmx_localtop_t *top;
2379 gmx_enerdata_t *enerd;
2380 gmx_global_stat_t gstat;
2381 t_graph *graph;
2382 real stepsize;
2383 real ustep;
2384 gmx_mdoutf_t outf;
2385 t_mdebin *mdebin;
2386 gmx_bool bDone, bAbort, do_x, do_f;
2387 tensor vir, pres;
2388 rvec mu_tot;
2389 int nsteps;
2390 int count = 0;
2391 int steps_accepted = 0;
2392 auto mdatoms = mdAtoms->mdatoms();
2394 /* Create 2 states on the stack and extract pointers that we will swap */
2395 em_state_t s0 {}, s1 {};
2396 em_state_t *s_min = &s0;
2397 em_state_t *s_try = &s1;
2399 /* Init em and store the local state in s_try */
2400 init_em(fplog, SD, cr, ms, outputProvider, inputrec, mdrunOptions,
2401 state_global, top_global, s_try, &top,
2402 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2403 vsite, constr, nullptr,
2404 nfile, fnm, &outf, &mdebin, wcycle);
2406 /* Print to log file */
2407 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2409 /* Set variables for stepsize (in nm). This is the largest
2410 * step that we are going to make in any direction.
2412 ustep = inputrec->em_stepsize;
2413 stepsize = 0;
2415 /* Max number of steps */
2416 nsteps = inputrec->nsteps;
2418 if (MASTER(cr))
2420 /* Print to the screen */
2421 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2423 if (fplog)
2425 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2427 EnergyEvaluator energyEvaluator {
2428 fplog, cr, ms,
2429 top_global, top,
2430 inputrec, nrnb, wcycle, gstat,
2431 vsite, constr, fcd, graph,
2432 mdAtoms, fr, enerd
2435 /**** HERE STARTS THE LOOP ****
2436 * count is the counter for the number of steps
2437 * bDone will be TRUE when the minimization has converged
2438 * bAbort will be TRUE when nsteps steps have been performed or when
2439 * the stepsize becomes smaller than is reasonable for machine precision
2441 count = 0;
2442 bDone = FALSE;
2443 bAbort = FALSE;
2444 while (!bDone && !bAbort)
2446 bAbort = (nsteps >= 0) && (count == nsteps);
2448 /* set new coordinates, except for first step */
2449 bool validStep = true;
2450 if (count > 0)
2452 validStep =
2453 do_em_step(cr, ms, inputrec, mdatoms, fr->bMolPBC,
2454 s_min, stepsize, &s_min->f, s_try,
2455 constr, top, nrnb, wcycle, count);
2458 if (validStep)
2460 energyEvaluator.run(s_try, mu_tot, vir, pres, count, count == 0);
2462 else
2464 // Signal constraint error during stepping with energy=inf
2465 s_try->epot = std::numeric_limits<real>::infinity();
2468 if (MASTER(cr))
2470 print_ebin_header(fplog, count, count);
2473 if (count == 0)
2475 s_min->epot = s_try->epot;
2478 /* Print it if necessary */
2479 if (MASTER(cr))
2481 if (mdrunOptions.verbose)
2483 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2484 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2485 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2486 fflush(stderr);
2489 if ( (count == 0) || (s_try->epot < s_min->epot) )
2491 /* Store the new (lower) energies */
2492 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2493 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2494 s_try->s.box, nullptr, nullptr, vir, pres, nullptr, mu_tot, constr);
2496 /* Prepare IMD energy record, if bIMD is TRUE. */
2497 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2499 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2500 do_per_step(steps_accepted, inputrec->nstdisreout),
2501 do_per_step(steps_accepted, inputrec->nstorireout),
2502 fplog, count, count, eprNORMAL,
2503 mdebin, fcd, &(top_global->groups), &(inputrec->opts), nullptr);
2504 fflush(fplog);
2508 /* Now if the new energy is smaller than the previous...
2509 * or if this is the first step!
2510 * or if we did random steps!
2513 if ( (count == 0) || (s_try->epot < s_min->epot) )
2515 steps_accepted++;
2517 /* Test whether the convergence criterion is met... */
2518 bDone = (s_try->fmax < inputrec->em_tol);
2520 /* Copy the arrays for force, positions and energy */
2521 /* The 'Min' array always holds the coords and forces of the minimal
2522 sampled energy */
2523 swap_em_state(&s_min, &s_try);
2524 if (count > 0)
2526 ustep *= 1.2;
2529 /* Write to trn, if necessary */
2530 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2531 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2532 write_em_traj(fplog, cr, outf, do_x, do_f, nullptr,
2533 top_global, inputrec, count,
2534 s_min, state_global, observablesHistory);
2536 else
2538 /* If energy is not smaller make the step smaller... */
2539 ustep *= 0.5;
2541 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2543 /* Reload the old state */
2544 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2545 s_min, top, mdAtoms, fr, vsite, constr,
2546 nrnb, wcycle);
2550 /* Determine new step */
2551 stepsize = ustep/s_min->fmax;
2553 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2554 #if GMX_DOUBLE
2555 if (count == nsteps || ustep < 1e-12)
2556 #else
2557 if (count == nsteps || ustep < 1e-6)
2558 #endif
2560 if (MASTER(cr))
2562 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != nullptr);
2563 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != nullptr);
2565 bAbort = TRUE;
2568 /* Send IMD energies and positions, if bIMD is TRUE. */
2569 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box,
2570 MASTER(cr) ? as_rvec_array(state_global->x.data()) : nullptr,
2571 inputrec, 0, wcycle) &&
2572 MASTER(cr))
2574 IMD_send_positions(inputrec->imd);
2577 count++;
2578 } /* End of the loop */
2580 /* IMD cleanup, if bIMD is TRUE. */
2581 IMD_finalize(inputrec->bIMD, inputrec->imd);
2583 /* Print some data... */
2584 if (MASTER(cr))
2586 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2588 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2589 top_global, inputrec, count,
2590 s_min, state_global, observablesHistory);
2592 if (MASTER(cr))
2594 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2596 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2597 s_min, sqrtNumAtoms);
2598 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2599 s_min, sqrtNumAtoms);
2602 finish_em(cr, outf, walltime_accounting, wcycle);
2604 /* To print the actual number of steps we needed somewhere */
2605 inputrec->nsteps = count;
2607 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2610 void
2611 Integrator::do_nm()
2613 const char *NM = "Normal Mode Analysis";
2614 gmx_mdoutf_t outf;
2615 int nnodes, node;
2616 gmx_localtop_t *top;
2617 gmx_enerdata_t *enerd;
2618 gmx_global_stat_t gstat;
2619 t_graph *graph;
2620 tensor vir, pres;
2621 rvec mu_tot;
2622 rvec *fneg, *dfdx;
2623 gmx_bool bSparse; /* use sparse matrix storage format */
2624 size_t sz;
2625 gmx_sparsematrix_t * sparse_matrix = nullptr;
2626 real * full_matrix = nullptr;
2628 /* added with respect to mdrun */
2629 int row, col;
2630 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2631 real x_min;
2632 bool bIsMaster = MASTER(cr);
2633 auto mdatoms = mdAtoms->mdatoms();
2635 if (constr != nullptr)
2637 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2640 gmx_shellfc_t *shellfc;
2642 em_state_t state_work {};
2644 /* Init em and store the local state in state_minimum */
2645 init_em(fplog, NM, cr, ms, outputProvider, inputrec, mdrunOptions,
2646 state_global, top_global, &state_work, &top,
2647 nrnb, mu_tot, fr, &enerd, &graph, mdAtoms, &gstat,
2648 vsite, constr, &shellfc,
2649 nfile, fnm, &outf, nullptr, wcycle);
2651 std::vector<size_t> atom_index = get_atom_index(top_global);
2652 snew(fneg, atom_index.size());
2653 snew(dfdx, atom_index.size());
2655 #if !GMX_DOUBLE
2656 if (bIsMaster)
2658 fprintf(stderr,
2659 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2660 " which MIGHT not be accurate enough for normal mode analysis.\n"
2661 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2662 " are fairly modest even if you recompile in double precision.\n\n");
2664 #endif
2666 /* Check if we can/should use sparse storage format.
2668 * Sparse format is only useful when the Hessian itself is sparse, which it
2669 * will be when we use a cutoff.
2670 * For small systems (n<1000) it is easier to always use full matrix format, though.
2672 if (EEL_FULL(fr->ic->eeltype) || fr->rlist == 0.0)
2674 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2675 bSparse = FALSE;
2677 else if (atom_index.size() < 1000)
2679 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%d), using full Hessian format.",
2680 atom_index.size());
2681 bSparse = FALSE;
2683 else
2685 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2686 bSparse = TRUE;
2689 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2690 sz = DIM*atom_index.size();
2692 fprintf(stderr, "Allocating Hessian memory...\n\n");
2694 if (bSparse)
2696 sparse_matrix = gmx_sparsematrix_init(sz);
2697 sparse_matrix->compressed_symmetric = TRUE;
2699 else
2701 snew(full_matrix, sz*sz);
2704 init_nrnb(nrnb);
2707 /* Write start time and temperature */
2708 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2710 /* fudge nr of steps to nr of atoms */
2711 inputrec->nsteps = atom_index.size()*2;
2713 if (bIsMaster)
2715 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2716 *(top_global->name), (int)inputrec->nsteps);
2719 nnodes = cr->nnodes;
2721 /* Make evaluate_energy do a single node force calculation */
2722 cr->nnodes = 1;
2723 EnergyEvaluator energyEvaluator {
2724 fplog, cr, ms,
2725 top_global, top,
2726 inputrec, nrnb, wcycle, gstat,
2727 vsite, constr, fcd, graph,
2728 mdAtoms, fr, enerd
2730 energyEvaluator.run(&state_work, mu_tot, vir, pres, -1, TRUE);
2731 cr->nnodes = nnodes;
2733 /* if forces are not small, warn user */
2734 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, &state_work);
2736 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work.fmax);
2737 if (state_work.fmax > 1.0e-3)
2739 GMX_LOG(mdlog.warning).appendText(
2740 "The force is probably not small enough to "
2741 "ensure that you are at a minimum.\n"
2742 "Be aware that negative eigenvalues may occur\n"
2743 "when the resulting matrix is diagonalized.");
2746 /***********************************************************
2748 * Loop over all pairs in matrix
2750 * do_force called twice. Once with positive and
2751 * once with negative displacement
2753 ************************************************************/
2755 /* Steps are divided one by one over the nodes */
2756 bool bNS = true;
2757 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2759 size_t atom = atom_index[aid];
2760 for (size_t d = 0; d < DIM; d++)
2762 gmx_int64_t step = 0;
2763 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2764 double t = 0;
2766 x_min = state_work.s.x[atom][d];
2768 for (unsigned int dx = 0; (dx < 2); dx++)
2770 if (dx == 0)
2772 state_work.s.x[atom][d] = x_min - der_range;
2774 else
2776 state_work.s.x[atom][d] = x_min + der_range;
2779 /* Make evaluate_energy do a single node force calculation */
2780 cr->nnodes = 1;
2781 if (shellfc)
2783 /* Now is the time to relax the shells */
2784 relax_shell_flexcon(fplog,
2787 mdrunOptions.verbose,
2788 step,
2789 inputrec,
2790 bNS,
2791 force_flags,
2792 top,
2793 constr,
2794 enerd,
2795 fcd,
2796 &state_work.s,
2797 state_work.f,
2798 vir,
2799 mdatoms,
2800 nrnb,
2801 wcycle,
2802 graph,
2803 &top_global->groups,
2804 shellfc,
2807 mu_tot,
2808 vsite,
2809 DdOpenBalanceRegionBeforeForceComputation::no,
2810 DdCloseBalanceRegionAfterForceComputation::no);
2811 bNS = false;
2812 step++;
2814 else
2816 energyEvaluator.run(&state_work, mu_tot, vir, pres, atom*2+dx, FALSE);
2819 cr->nnodes = nnodes;
2821 if (dx == 0)
2823 for (size_t i = 0; i < atom_index.size(); i++)
2825 copy_rvec(state_work.f[atom_index[i]], fneg[i]);
2830 /* x is restored to original */
2831 state_work.s.x[atom][d] = x_min;
2833 for (size_t j = 0; j < atom_index.size(); j++)
2835 for (size_t k = 0; (k < DIM); k++)
2837 dfdx[j][k] =
2838 -(state_work.f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
2842 if (!bIsMaster)
2844 #if GMX_MPI
2845 #define mpi_type GMX_MPI_REAL
2846 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
2847 cr->nodeid, cr->mpi_comm_mygroup);
2848 #endif
2850 else
2852 for (node = 0; (node < nnodes && atom+node < atom_index.size()); node++)
2854 if (node > 0)
2856 #if GMX_MPI
2857 MPI_Status stat;
2858 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
2859 cr->mpi_comm_mygroup, &stat);
2860 #undef mpi_type
2861 #endif
2864 row = (atom + node)*DIM + d;
2866 for (size_t j = 0; j < atom_index.size(); j++)
2868 for (size_t k = 0; k < DIM; k++)
2870 col = j*DIM + k;
2872 if (bSparse)
2874 if (col >= row && dfdx[j][k] != 0.0)
2876 gmx_sparsematrix_increment_value(sparse_matrix,
2877 row, col, dfdx[j][k]);
2880 else
2882 full_matrix[row*sz+col] = dfdx[j][k];
2889 if (mdrunOptions.verbose && fplog)
2891 fflush(fplog);
2894 /* write progress */
2895 if (bIsMaster && mdrunOptions.verbose)
2897 fprintf(stderr, "\rFinished step %d out of %d",
2898 static_cast<int>(std::min(atom+nnodes, atom_index.size())),
2899 static_cast<int>(atom_index.size()));
2900 fflush(stderr);
2904 if (bIsMaster)
2906 fprintf(stderr, "\n\nWriting Hessian...\n");
2907 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2910 finish_em(cr, outf, walltime_accounting, wcycle);
2912 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);
2915 } // namespace gmx