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, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file defines integrators for energy minimization
41 * \author Berk Hess <hess@kth.se>
42 * \author Erik Lindahl <erik@kth.se>
43 * \ingroup module_mdlib
58 #include "gromacs/commandline/filenm.h"
59 #include "gromacs/domdec/domdec.h"
60 #include "gromacs/domdec/domdec_struct.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/fileio/confio.h"
63 #include "gromacs/fileio/mtxio.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/imd/imd.h"
67 #include "gromacs/linearalgebra/sparsematrix.h"
68 #include "gromacs/listed-forces/manage-threading.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/vec.h"
71 #include "gromacs/mdlib/constr.h"
72 #include "gromacs/mdlib/force.h"
73 #include "gromacs/mdlib/forcerec.h"
74 #include "gromacs/mdlib/gmx_omp_nthreads.h"
75 #include "gromacs/mdlib/md_support.h"
76 #include "gromacs/mdlib/mdatoms.h"
77 #include "gromacs/mdlib/mdebin.h"
78 #include "gromacs/mdlib/mdrun.h"
79 #include "gromacs/mdlib/mdsetup.h"
80 #include "gromacs/mdlib/ns.h"
81 #include "gromacs/mdlib/shellfc.h"
82 #include "gromacs/mdlib/sim_util.h"
83 #include "gromacs/mdlib/tgroup.h"
84 #include "gromacs/mdlib/trajectory_writing.h"
85 #include "gromacs/mdlib/update.h"
86 #include "gromacs/mdlib/vsite.h"
87 #include "gromacs/mdtypes/commrec.h"
88 #include "gromacs/mdtypes/inputrec.h"
89 #include "gromacs/mdtypes/md_enums.h"
90 #include "gromacs/pbcutil/mshift.h"
91 #include "gromacs/pbcutil/pbc.h"
92 #include "gromacs/timing/wallcycle.h"
93 #include "gromacs/timing/walltime_accounting.h"
94 #include "gromacs/topology/mtop_util.h"
95 #include "gromacs/topology/topology.h"
96 #include "gromacs/utility/cstringutil.h"
97 #include "gromacs/utility/exceptions.h"
98 #include "gromacs/utility/fatalerror.h"
99 #include "gromacs/utility/logger.h"
100 #include "gromacs/utility/smalloc.h"
102 //! Utility structure for manipulating states during EM
104 //! Copy of the global state
110 //! Norm of the force
118 //! Initiate em_state_t structure and return pointer to it
119 static em_state_t *init_em_state()
125 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
126 snew(ems->s.lambda, efptNR);
131 //! Print the EM starting conditions
132 static void print_em_start(FILE *fplog,
134 gmx_walltime_accounting_t walltime_accounting,
135 gmx_wallcycle_t wcycle,
138 walltime_accounting_start(walltime_accounting);
139 wallcycle_start(wcycle, ewcRUN);
140 print_start(fplog, cr, walltime_accounting, name);
143 //! Stop counting time for EM
144 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
145 gmx_wallcycle_t wcycle)
147 wallcycle_stop(wcycle, ewcRUN);
149 walltime_accounting_end(walltime_accounting);
152 //! Printing a log file and console header
153 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
156 fprintf(out, "%s:\n", minimizer);
157 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
158 fprintf(out, " Number of steps = %12d\n", nsteps);
161 //! Print warning message
162 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
168 "\nEnergy minimization reached the maximum number "
169 "of steps before the forces reached the requested "
170 "precision Fmax < %g.\n", ftol);
175 "\nEnergy minimization has stopped, but the forces have "
176 "not converged to the requested precision Fmax < %g (which "
177 "may not be possible for your system). It stopped "
178 "because the algorithm tried to make a new step whose size "
179 "was too small, or there was no change in the energy since "
180 "last step. Either way, we regard the minimization as "
181 "converged to within the available machine precision, "
182 "given your starting configuration and EM parameters.\n%s%s",
184 sizeof(real) < sizeof(double) ?
185 "\nDouble precision normally gives you higher accuracy, but "
186 "this is often not needed for preparing to run molecular "
190 "You might need to increase your constraint accuracy, or turn\n"
191 "off constraints altogether (set constraints = none in mdp file)\n" :
194 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
197 //! Print message about convergence of the EM
198 static void print_converged(FILE *fp, const char *alg, real ftol,
199 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
200 real epot, real fmax, int nfmax, real fnorm)
202 char buf[STEPSTRSIZE];
206 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
207 alg, ftol, gmx_step_str(count, buf));
209 else if (count < nsteps)
211 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
212 "but did not reach the requested Fmax < %g.\n",
213 alg, gmx_step_str(count, buf), ftol);
217 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
218 alg, ftol, gmx_step_str(count, buf));
222 fprintf(fp, "Potential Energy = %21.14e\n", epot);
223 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
224 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
226 fprintf(fp, "Potential Energy = %14.7e\n", epot);
227 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
228 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
232 //! Compute the norm and max of the force array in parallel
233 static void get_f_norm_max(t_commrec *cr,
234 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
235 real *fnorm, real *fmax, int *a_fmax)
239 int la_max, a_max, start, end, i, m, gf;
241 /* This routine finds the largest force and returns it.
242 * On parallel machines the global max is taken.
248 end = mdatoms->homenr;
249 if (mdatoms->cFREEZE)
251 for (i = start; i < end; i++)
253 gf = mdatoms->cFREEZE[i];
255 for (m = 0; m < DIM; m++)
257 if (!opts->nFreeze[gf][m])
259 fam += gmx::square(f[i][m]);
272 for (i = start; i < end; i++)
284 if (la_max >= 0 && DOMAINDECOMP(cr))
286 a_max = cr->dd->gatindex[la_max];
294 snew(sum, 2*cr->nnodes+1);
295 sum[2*cr->nodeid] = fmax2;
296 sum[2*cr->nodeid+1] = a_max;
297 sum[2*cr->nnodes] = fnorm2;
298 gmx_sumd(2*cr->nnodes+1, sum, cr);
299 fnorm2 = sum[2*cr->nnodes];
300 /* Determine the global maximum */
301 for (i = 0; i < cr->nnodes; i++)
303 if (sum[2*i] > fmax2)
306 a_max = (int)(sum[2*i+1] + 0.5);
314 *fnorm = sqrt(fnorm2);
326 //! Compute the norm of the force
327 static void get_state_f_norm_max(t_commrec *cr,
328 t_grpopts *opts, t_mdatoms *mdatoms,
331 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
334 //! Initialize the energy minimization
335 void init_em(FILE *fplog, const char *title,
336 t_commrec *cr, t_inputrec *ir,
337 t_state *state_global, gmx_mtop_t *top_global,
338 em_state_t *ems, gmx_localtop_t **top,
340 t_nrnb *nrnb, rvec mu_tot,
341 t_forcerec *fr, gmx_enerdata_t **enerd,
342 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
343 gmx_vsite_t *vsite, gmx_constr_t constr, gmx_shellfc_t **shellfc,
344 int nfile, const t_filenm fnm[],
345 gmx_mdoutf_t *outf, t_mdebin **mdebin,
346 int imdport, unsigned long gmx_unused Flags,
347 gmx_wallcycle_t wcycle)
354 fprintf(fplog, "Initiating %s\n", title);
357 state_global->ngtc = 0;
359 /* Initialize lambda variables */
360 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
364 /* Interactive molecular dynamics */
365 init_IMD(ir, cr, top_global, fplog, 1, state_global->x,
366 nfile, fnm, NULL, imdport, Flags);
370 GMX_ASSERT(shellfc != NULL, "With NM we always support shells");
372 *shellfc = init_shell_flexcon(stdout,
374 n_flexible_constraints(constr),
380 GMX_ASSERT(EI_ENERGY_MINIMIZATION(ir->eI), "This else currently only handles energy minimizers, consider if your algorithm needs shell/flexible-constraint support");
382 /* With energy minimization, shells and flexible constraints are
383 * automatically minimized when treated like normal DOFS.
391 if (DOMAINDECOMP(cr))
393 *top = dd_init_local_top(top_global);
395 dd_init_local_state(cr->dd, state_global, &ems->s);
399 /* Distribute the charge groups over the nodes from the master node */
400 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
401 state_global, top_global, ir,
402 &ems->s, &ems->f, mdatoms, *top,
405 dd_store_state(cr->dd, &ems->s);
411 snew(*f, top_global->natoms);
413 /* Just copy the state */
414 ems->s = *state_global;
415 /* We need to allocate one element extra, since we might use
416 * (unaligned) 4-wide SIMD loads to access rvec entries.
418 snew(ems->s.x, ems->s.nalloc + 1);
419 snew(ems->f, ems->s.nalloc+1);
420 snew(ems->s.v, ems->s.nalloc+1);
421 for (i = 0; i < state_global->natoms; i++)
423 copy_rvec(state_global->x[i], ems->s.x[i]);
425 copy_mat(state_global->box, ems->s.box);
429 mdAlgorithmsSetupAtomData(cr, ir, top_global, *top, fr,
431 vsite, shellfc ? *shellfc : NULL);
433 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
437 set_vsite_top(vsite, *top, mdatoms, cr);
443 if (ir->eConstrAlg == econtSHAKE &&
444 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
446 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
447 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
450 if (!DOMAINDECOMP(cr))
452 set_constraints(constr, *top, ir, mdatoms, cr);
455 if (!ir->bContinuation)
457 /* Constrain the starting coordinates */
459 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
460 ir, cr, -1, 0, 1.0, mdatoms,
461 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
462 ems->s.lambda[efptFEP], &dvdl_constr,
463 NULL, NULL, nrnb, econqCoord);
469 *gstat = global_stat_init(ir);
476 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, NULL, wcycle);
479 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
484 /* Init bin for energy stuff */
485 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
489 calc_shifts(ems->s.box, fr->shift_vec);
492 //! Finalize the minimization
493 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
494 gmx_walltime_accounting_t walltime_accounting,
495 gmx_wallcycle_t wcycle)
497 if (!(cr->duty & DUTY_PME))
499 /* Tell the PME only node to finish */
500 gmx_pme_send_finish(cr);
505 em_time_end(walltime_accounting, wcycle);
508 //! Swap two different EM states during minimization
509 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
518 //! Copy coordinate from an EM state to a "normal" state structure
519 static void copy_em_coords(em_state_t *ems, t_state *state)
523 for (i = 0; (i < state->natoms); i++)
525 copy_rvec(ems->s.x[i], state->x[i]);
529 //! Save the EM trajectory
530 static void write_em_traj(FILE *fplog, t_commrec *cr,
532 gmx_bool bX, gmx_bool bF, const char *confout,
533 gmx_mtop_t *top_global,
534 t_inputrec *ir, gmx_int64_t step,
536 t_state *state_global)
539 gmx_bool bIMDout = FALSE;
542 /* Shall we do IMD output? */
545 bIMDout = do_per_step(step, IMD_get_step(ir->imd->setup));
548 if ((bX || bF || bIMDout || confout != NULL) && !DOMAINDECOMP(cr))
550 copy_em_coords(state, state_global);
556 mdof_flags |= MDOF_X;
560 mdof_flags |= MDOF_F;
563 /* If we want IMD output, set appropriate MDOF flag */
566 mdof_flags |= MDOF_IMD;
569 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
570 top_global, step, (double)step,
571 &state->s, state_global, state->f);
573 if (confout != NULL && MASTER(cr))
575 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
577 /* Make molecules whole only for confout writing */
578 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
582 write_sto_conf_mtop(confout,
583 *top_global->name, top_global,
584 state_global->x, NULL, ir->ePBC, state_global->box);
588 //! \brief Do one minimization step
590 // \returns true when the step succeeded, false when a constraint error occurred
591 static bool do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
593 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
594 gmx_constr_t constr, gmx_localtop_t *top,
595 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
602 int nthreads gmx_unused;
604 bool validStep = true;
609 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
611 gmx_incons("state mismatch in do_em_step");
614 s2->flags = s1->flags;
616 if (s2->nalloc != s1->nalloc)
618 s2->nalloc = s1->nalloc;
619 /* We need to allocate one element extra, since we might use
620 * (unaligned) 4-wide SIMD loads to access rvec entries.
622 srenew(s2->x, s1->nalloc + 1);
623 srenew(ems2->f, s1->nalloc);
624 if (s2->flags & (1<<estCGP))
626 srenew(s2->cg_p, s1->nalloc + 1);
630 s2->natoms = s1->natoms;
631 copy_mat(s1->box, s2->box);
632 /* Copy free energy state */
633 for (int i = 0; i < efptNR; i++)
635 s2->lambda[i] = s1->lambda[i];
637 copy_mat(s1->box, s2->box);
642 // cppcheck-suppress unreadVariable
643 nthreads = gmx_omp_nthreads_get(emntUpdate);
644 #pragma omp parallel num_threads(nthreads)
650 #pragma omp for schedule(static) nowait
651 for (int i = start; i < end; i++)
659 for (int m = 0; m < DIM; m++)
661 if (ir->opts.nFreeze[gf][m])
667 x2[i][m] = x1[i][m] + a*f[i][m];
671 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
674 if (s2->flags & (1<<estCGP))
676 /* Copy the CG p vector */
679 #pragma omp for schedule(static) nowait
680 for (int i = start; i < end; i++)
682 // Trivial OpenMP block that does not throw
683 copy_rvec(p1[i], p2[i]);
687 if (DOMAINDECOMP(cr))
689 s2->ddp_count = s1->ddp_count;
690 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
693 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
696 /* We need to allocate one element extra, since we might use
697 * (unaligned) 4-wide SIMD loads to access rvec entries.
699 srenew(s2->cg_gl, s2->cg_gl_nalloc + 1);
701 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
704 s2->ncg_gl = s1->ncg_gl;
705 #pragma omp for schedule(static) nowait
706 for (int i = 0; i < s2->ncg_gl; i++)
708 s2->cg_gl[i] = s1->cg_gl[i];
710 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
716 wallcycle_start(wcycle, ewcCONSTR);
719 constrain(NULL, TRUE, TRUE, constr, &top->idef,
720 ir, cr, count, 0, 1.0, md,
721 s1->x, s2->x, NULL, bMolPBC, s2->box,
722 s2->lambda[efptBONDED], &dvdl_constr,
723 NULL, NULL, nrnb, econqCoord);
724 wallcycle_stop(wcycle, ewcCONSTR);
726 // We should move this check to the different minimizers
727 if (!validStep && ir->eI != eiSteep)
729 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
730 EI(ir->eI), EI(eiSteep), EI(ir->eI));
737 //! Prepare EM for using domain decomposition parallellization
738 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
739 gmx_mtop_t *top_global, t_inputrec *ir,
740 em_state_t *ems, gmx_localtop_t *top,
741 t_mdatoms *mdatoms, t_forcerec *fr,
742 gmx_vsite_t *vsite, gmx_constr_t constr,
743 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
745 /* Repartition the domain decomposition */
746 dd_partition_system(fplog, step, cr, FALSE, 1,
747 NULL, top_global, ir,
749 mdatoms, top, fr, vsite, constr,
750 nrnb, wcycle, FALSE);
751 dd_store_state(cr->dd, &ems->s);
754 //! De one energy evaluation
755 static void evaluate_energy(FILE *fplog, t_commrec *cr,
756 gmx_mtop_t *top_global,
757 em_state_t *ems, gmx_localtop_t *top,
758 t_inputrec *inputrec,
759 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
760 gmx_global_stat_t gstat,
761 gmx_vsite_t *vsite, gmx_constr_t constr,
763 t_graph *graph, t_mdatoms *mdatoms,
764 t_forcerec *fr, rvec mu_tot,
765 gmx_enerdata_t *enerd, tensor vir, tensor pres,
766 gmx_int64_t count, gmx_bool bFirst)
770 tensor force_vir, shake_vir, ekin;
771 real dvdl_constr, prescorr, enercorr, dvdlcorr;
774 /* Set the time to the initial time, the time does not change during EM */
775 t = inputrec->init_t;
778 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
780 /* This is the first state or an old state used before the last ns */
786 if (inputrec->nstlist > 0)
794 construct_vsites(vsite, ems->s.x, 1, NULL,
795 top->idef.iparams, top->idef.il,
796 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
799 if (DOMAINDECOMP(cr) && bNS)
801 /* Repartition the domain decomposition */
802 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
803 ems, top, mdatoms, fr, vsite, constr,
807 /* Calc force & energy on new trial position */
808 /* do_force always puts the charge groups in the box and shifts again
809 * We do not unshift, so molecules are always whole in congrad.c
811 do_force(fplog, cr, inputrec,
812 count, nrnb, wcycle, top, &top_global->groups,
813 ems->s.box, ems->s.x, &ems->s.hist,
814 ems->f, force_vir, mdatoms, enerd, fcd,
815 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
816 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
817 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
818 (bNS ? GMX_FORCE_NS : 0));
820 /* Clear the unused shake virial and pressure */
821 clear_mat(shake_vir);
824 /* Communicate stuff when parallel */
825 if (PAR(cr) && inputrec->eI != eiNM)
827 wallcycle_start(wcycle, ewcMoveE);
829 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
830 inputrec, NULL, NULL, NULL, 1, &terminate,
836 wallcycle_stop(wcycle, ewcMoveE);
839 /* Calculate long range corrections to pressure and energy */
840 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
841 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
842 enerd->term[F_DISPCORR] = enercorr;
843 enerd->term[F_EPOT] += enercorr;
844 enerd->term[F_PRES] += prescorr;
845 enerd->term[F_DVDL] += dvdlcorr;
847 ems->epot = enerd->term[F_EPOT];
851 /* Project out the constraint components of the force */
852 wallcycle_start(wcycle, ewcCONSTR);
854 constrain(NULL, FALSE, FALSE, constr, &top->idef,
855 inputrec, cr, count, 0, 1.0, mdatoms,
856 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
857 ems->s.lambda[efptBONDED], &dvdl_constr,
858 NULL, &shake_vir, nrnb, econqForceDispl);
859 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
860 m_add(force_vir, shake_vir, vir);
861 wallcycle_stop(wcycle, ewcCONSTR);
865 copy_mat(force_vir, vir);
869 enerd->term[F_PRES] =
870 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
872 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
874 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
876 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
880 //! Parallel utility summing energies and forces
881 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
882 gmx_mtop_t *top_global,
883 em_state_t *s_min, em_state_t *s_b)
887 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
889 unsigned char *grpnrFREEZE;
893 fprintf(debug, "Doing reorder_partsum\n");
899 cgs_gl = dd_charge_groups_global(cr->dd);
900 index = cgs_gl->index;
902 /* Collect fm in a global vector fmg.
903 * This conflicts with the spirit of domain decomposition,
904 * but to fully optimize this a much more complicated algorithm is required.
906 snew(fmg, top_global->natoms);
908 ncg = s_min->s.ncg_gl;
909 cg_gl = s_min->s.cg_gl;
911 for (c = 0; c < ncg; c++)
916 for (a = a0; a < a1; a++)
918 copy_rvec(fm[i], fmg[a]);
922 gmx_sum(top_global->natoms*3, fmg[0], cr);
924 /* Now we will determine the part of the sum for the cgs in state s_b */
926 cg_gl = s_b->s.cg_gl;
930 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
931 for (c = 0; c < ncg; c++)
936 for (a = a0; a < a1; a++)
938 if (mdatoms->cFREEZE && grpnrFREEZE)
942 for (m = 0; m < DIM; m++)
944 if (!opts->nFreeze[gf][m])
946 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
958 //! Print some stuff, like beta, whatever that means.
959 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
960 gmx_mtop_t *top_global,
961 em_state_t *s_min, em_state_t *s_b)
967 /* This is just the classical Polak-Ribiere calculation of beta;
968 * it looks a bit complicated since we take freeze groups into account,
969 * and might have to sum it in parallel runs.
972 if (!DOMAINDECOMP(cr) ||
973 (s_min->s.ddp_count == cr->dd->ddp_count &&
974 s_b->s.ddp_count == cr->dd->ddp_count))
980 /* This part of code can be incorrect with DD,
981 * since the atom ordering in s_b and s_min might differ.
983 for (i = 0; i < mdatoms->homenr; i++)
985 if (mdatoms->cFREEZE)
987 gf = mdatoms->cFREEZE[i];
989 for (m = 0; m < DIM; m++)
991 if (!opts->nFreeze[gf][m])
993 sum += (fb[i][m] - fm[i][m])*fb[i][m];
1000 /* We need to reorder cgs while summing */
1001 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
1005 gmx_sumd(1, &sum, cr);
1008 return sum/gmx::square(s_min->fnorm);
1014 /*! \brief Do conjugate gradients minimization
1015 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
1016 int nfile, const t_filenm fnm[],
1017 const gmx_output_env_t *oenv, gmx_bool bVerbose,
1019 gmx_vsite_t *vsite, gmx_constr_t constr,
1021 t_inputrec *inputrec,
1022 gmx_mtop_t *top_global, t_fcdata *fcd,
1023 t_state *state_global,
1025 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1028 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1029 gmx_membed_t gmx_unused *membed,
1030 real cpt_period, real max_hours,
1032 unsigned long Flags,
1033 gmx_walltime_accounting_t walltime_accounting)
1035 double do_cg(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
1036 int nfile, const t_filenm fnm[],
1037 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
1038 int gmx_unused nstglobalcomm,
1039 gmx_vsite_t *vsite, gmx_constr_t constr,
1040 int gmx_unused stepout,
1041 t_inputrec *inputrec,
1042 gmx_mtop_t *top_global, t_fcdata *fcd,
1043 t_state *state_global,
1045 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1046 gmx_edsam_t gmx_unused ed,
1048 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1049 gmx_membed_t gmx_unused *membed,
1050 real gmx_unused cpt_period, real gmx_unused max_hours,
1052 unsigned long gmx_unused Flags,
1053 gmx_walltime_accounting_t walltime_accounting)
1055 const char *CG = "Polak-Ribiere Conjugate Gradients";
1057 em_state_t *s_min, *s_a, *s_b, *s_c;
1058 gmx_localtop_t *top;
1059 gmx_enerdata_t *enerd;
1061 gmx_global_stat_t gstat;
1064 double gpa, gpb, gpc, tmp, minstep;
1067 real a, b, c, beta = 0.0;
1071 gmx_bool converged, foundlower;
1073 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1075 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1077 int i, m, gf, step, nminstep;
1081 s_min = init_em_state();
1082 s_a = init_em_state();
1083 s_b = init_em_state();
1084 s_c = init_em_state();
1086 /* Init em and store the local state in s_min */
1087 init_em(fplog, CG, cr, inputrec,
1088 state_global, top_global, s_min, &top, &f,
1089 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
1090 vsite, constr, NULL,
1091 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1093 /* Print to log file */
1094 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1096 /* Max number of steps */
1097 number_steps = inputrec->nsteps;
1101 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1105 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1108 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1109 /* do_force always puts the charge groups in the box and shifts again
1110 * We do not unshift, so molecules are always whole in congrad.c
1112 evaluate_energy(fplog, cr,
1113 top_global, s_min, top,
1114 inputrec, nrnb, wcycle, gstat,
1115 vsite, constr, fcd, graph, mdatoms, fr,
1116 mu_tot, enerd, vir, pres, -1, TRUE);
1121 /* Copy stuff to the energy bin for easy printing etc. */
1122 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1123 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1124 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1126 print_ebin_header(fplog, step, step);
1127 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1128 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1132 /* Estimate/guess the initial stepsize */
1133 stepsize = inputrec->em_stepsize/s_min->fnorm;
1137 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1138 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1139 s_min->fmax, s_min->a_fmax+1);
1140 fprintf(stderr, " F-Norm = %12.5e\n",
1141 s_min->fnorm/sqrtNumAtoms);
1142 fprintf(stderr, "\n");
1143 /* and copy to the log file too... */
1144 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1145 s_min->fmax, s_min->a_fmax+1);
1146 fprintf(fplog, " F-Norm = %12.5e\n",
1147 s_min->fnorm/sqrtNumAtoms);
1148 fprintf(fplog, "\n");
1150 /* Start the loop over CG steps.
1151 * Each successful step is counted, and we continue until
1152 * we either converge or reach the max number of steps.
1155 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1158 /* start taking steps in a new direction
1159 * First time we enter the routine, beta=0, and the direction is
1160 * simply the negative gradient.
1163 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1168 for (i = 0; i < mdatoms->homenr; i++)
1170 if (mdatoms->cFREEZE)
1172 gf = mdatoms->cFREEZE[i];
1174 for (m = 0; m < DIM; m++)
1176 if (!inputrec->opts.nFreeze[gf][m])
1178 p[i][m] = sf[i][m] + beta*p[i][m];
1179 gpa -= p[i][m]*sf[i][m];
1180 /* f is negative gradient, thus the sign */
1189 /* Sum the gradient along the line across CPUs */
1192 gmx_sumd(1, &gpa, cr);
1195 /* Calculate the norm of the search vector */
1196 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1198 /* Just in case stepsize reaches zero due to numerical precision... */
1201 stepsize = inputrec->em_stepsize/pnorm;
1205 * Double check the value of the derivative in the search direction.
1206 * If it is positive it must be due to the old information in the
1207 * CG formula, so just remove that and start over with beta=0.
1208 * This corresponds to a steepest descent step.
1213 step--; /* Don't count this step since we are restarting */
1214 continue; /* Go back to the beginning of the big for-loop */
1217 /* Calculate minimum allowed stepsize, before the average (norm)
1218 * relative change in coordinate is smaller than precision
1221 for (i = 0; i < mdatoms->homenr; i++)
1223 for (m = 0; m < DIM; m++)
1225 tmp = fabs(s_min->s.x[i][m]);
1234 /* Add up from all CPUs */
1237 gmx_sumd(1, &minstep, cr);
1240 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1242 if (stepsize < minstep)
1248 /* Write coordinates if necessary */
1249 do_x = do_per_step(step, inputrec->nstxout);
1250 do_f = do_per_step(step, inputrec->nstfout);
1252 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1253 top_global, inputrec, step,
1254 s_min, state_global);
1256 /* Take a step downhill.
1257 * In theory, we should minimize the function along this direction.
1258 * That is quite possible, but it turns out to take 5-10 function evaluations
1259 * for each line. However, we dont really need to find the exact minimum -
1260 * it is much better to start a new CG step in a modified direction as soon
1261 * as we are close to it. This will save a lot of energy evaluations.
1263 * In practice, we just try to take a single step.
1264 * If it worked (i.e. lowered the energy), we increase the stepsize but
1265 * the continue straight to the next CG step without trying to find any minimum.
1266 * If it didn't work (higher energy), there must be a minimum somewhere between
1267 * the old position and the new one.
1269 * Due to the finite numerical accuracy, it turns out that it is a good idea
1270 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1271 * This leads to lower final energies in the tests I've done. / Erik
1273 s_a->epot = s_min->epot;
1275 c = a + stepsize; /* reference position along line is zero */
1277 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1279 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1280 s_min, top, mdatoms, fr, vsite, constr,
1284 /* Take a trial step (new coords in s_c) */
1285 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1286 constr, top, nrnb, wcycle, -1);
1289 /* Calculate energy for the trial step */
1290 evaluate_energy(fplog, cr,
1291 top_global, s_c, top,
1292 inputrec, nrnb, wcycle, gstat,
1293 vsite, constr, fcd, graph, mdatoms, fr,
1294 mu_tot, enerd, vir, pres, -1, FALSE);
1296 /* Calc derivative along line */
1300 for (i = 0; i < mdatoms->homenr; i++)
1302 for (m = 0; m < DIM; m++)
1304 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1307 /* Sum the gradient along the line across CPUs */
1310 gmx_sumd(1, &gpc, cr);
1313 /* This is the max amount of increase in energy we tolerate */
1314 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1316 /* Accept the step if the energy is lower, or if it is not significantly higher
1317 * and the line derivative is still negative.
1319 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1322 /* Great, we found a better energy. Increase step for next iteration
1323 * if we are still going down, decrease it otherwise
1327 stepsize *= 1.618034; /* The golden section */
1331 stepsize *= 0.618034; /* 1/golden section */
1336 /* New energy is the same or higher. We will have to do some work
1337 * to find a smaller value in the interval. Take smaller step next time!
1340 stepsize *= 0.618034;
1346 /* OK, if we didn't find a lower value we will have to locate one now - there must
1347 * be one in the interval [a=0,c].
1348 * The same thing is valid here, though: Don't spend dozens of iterations to find
1349 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1350 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1352 * I also have a safeguard for potentially really pathological functions so we never
1353 * take more than 20 steps before we give up ...
1355 * If we already found a lower value we just skip this step and continue to the update.
1363 /* Select a new trial point.
1364 * If the derivatives at points a & c have different sign we interpolate to zero,
1365 * otherwise just do a bisection.
1367 if (gpa < 0 && gpc > 0)
1369 b = a + gpa*(a-c)/(gpc-gpa);
1376 /* safeguard if interpolation close to machine accuracy causes errors:
1377 * never go outside the interval
1379 if (b <= a || b >= c)
1384 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1386 /* Reload the old state */
1387 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1388 s_min, top, mdatoms, fr, vsite, constr,
1392 /* Take a trial step to this new point - new coords in s_b */
1393 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1394 constr, top, nrnb, wcycle, -1);
1397 /* Calculate energy for the trial step */
1398 evaluate_energy(fplog, cr,
1399 top_global, s_b, top,
1400 inputrec, nrnb, wcycle, gstat,
1401 vsite, constr, fcd, graph, mdatoms, fr,
1402 mu_tot, enerd, vir, pres, -1, FALSE);
1404 /* p does not change within a step, but since the domain decomposition
1405 * might change, we have to use cg_p of s_b here.
1410 for (i = 0; i < mdatoms->homenr; i++)
1412 for (m = 0; m < DIM; m++)
1414 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1417 /* Sum the gradient along the line across CPUs */
1420 gmx_sumd(1, &gpb, cr);
1425 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1426 s_a->epot, s_b->epot, s_c->epot, gpb);
1429 epot_repl = s_b->epot;
1431 /* Keep one of the intervals based on the value of the derivative at the new point */
1434 /* Replace c endpoint with b */
1435 swap_em_state(s_b, s_c);
1441 /* Replace a endpoint with b */
1442 swap_em_state(s_b, s_a);
1448 * Stop search as soon as we find a value smaller than the endpoints.
1449 * Never run more than 20 steps, no matter what.
1453 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1456 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1459 /* OK. We couldn't find a significantly lower energy.
1460 * If beta==0 this was steepest descent, and then we give up.
1461 * If not, set beta=0 and restart with steepest descent before quitting.
1471 /* Reset memory before giving up */
1477 /* Select min energy state of A & C, put the best in B.
1479 if (s_c->epot < s_a->epot)
1483 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1484 s_c->epot, s_a->epot);
1486 swap_em_state(s_b, s_c);
1493 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1494 s_a->epot, s_c->epot);
1496 swap_em_state(s_b, s_a);
1505 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1508 swap_em_state(s_b, s_c);
1512 /* new search direction */
1513 /* beta = 0 means forget all memory and restart with steepest descents. */
1514 if (nstcg && ((step % nstcg) == 0))
1520 /* s_min->fnorm cannot be zero, because then we would have converged
1524 /* Polak-Ribiere update.
1525 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1527 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1529 /* Limit beta to prevent oscillations */
1530 if (fabs(beta) > 5.0)
1536 /* update positions */
1537 swap_em_state(s_min, s_b);
1540 /* Print it if necessary */
1545 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1546 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1547 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1548 s_min->fmax, s_min->a_fmax+1);
1551 /* Store the new (lower) energies */
1552 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1553 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1554 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1556 do_log = do_per_step(step, inputrec->nstlog);
1557 do_ene = do_per_step(step, inputrec->nstenergy);
1559 /* Prepare IMD energy record, if bIMD is TRUE. */
1560 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1564 print_ebin_header(fplog, step, step);
1566 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1567 do_log ? fplog : NULL, step, step, eprNORMAL,
1568 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1571 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1572 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1574 IMD_send_positions(inputrec->imd);
1577 /* Stop when the maximum force lies below tolerance.
1578 * If we have reached machine precision, converged is already set to true.
1580 converged = converged || (s_min->fmax < inputrec->em_tol);
1582 } /* End of the loop */
1584 /* IMD cleanup, if bIMD is TRUE. */
1585 IMD_finalize(inputrec->bIMD, inputrec->imd);
1589 step--; /* we never took that last step in this case */
1592 if (s_min->fmax > inputrec->em_tol)
1596 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1597 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1604 /* If we printed energy and/or logfile last step (which was the last step)
1605 * we don't have to do it again, but otherwise print the final values.
1609 /* Write final value to log since we didn't do anything the last step */
1610 print_ebin_header(fplog, step, step);
1612 if (!do_ene || !do_log)
1614 /* Write final energy file entries */
1615 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1616 !do_log ? fplog : NULL, step, step, eprNORMAL,
1617 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1621 /* Print some stuff... */
1624 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1628 * For accurate normal mode calculation it is imperative that we
1629 * store the last conformation into the full precision binary trajectory.
1631 * However, we should only do it if we did NOT already write this step
1632 * above (which we did if do_x or do_f was true).
1634 do_x = !do_per_step(step, inputrec->nstxout);
1635 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1637 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1638 top_global, inputrec, step,
1639 s_min, state_global);
1644 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1645 fnormn = s_min->fnorm/sqrtNumAtoms;
1646 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1647 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1648 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1649 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1651 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1654 finish_em(cr, outf, walltime_accounting, wcycle);
1656 /* To print the actual number of steps we needed somewhere */
1657 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1660 } /* That's all folks */
1663 /*! \brief Do L-BFGS conjugate gradients minimization
1664 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
1665 int nfile, const t_filenm fnm[],
1666 const gmx_output_env_t *oenv, gmx_bool bVerbose,
1668 gmx_vsite_t *vsite, gmx_constr_t constr,
1670 t_inputrec *inputrec,
1671 gmx_mtop_t *top_global, t_fcdata *fcd,
1672 t_state *state_global,
1674 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1677 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1678 real cpt_period, real max_hours,
1680 unsigned long Flags,
1681 gmx_walltime_accounting_t walltime_accounting)
1683 double do_lbfgs(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
1684 int nfile, const t_filenm fnm[],
1685 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
1686 int gmx_unused nstglobalcomm,
1687 gmx_vsite_t *vsite, gmx_constr_t constr,
1688 int gmx_unused stepout,
1689 t_inputrec *inputrec,
1690 gmx_mtop_t *top_global, t_fcdata *fcd,
1691 t_state *state_global,
1693 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1694 gmx_edsam_t gmx_unused ed,
1696 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1697 gmx_membed_t gmx_unused *membed,
1698 real gmx_unused cpt_period, real gmx_unused max_hours,
1700 unsigned long gmx_unused Flags,
1701 gmx_walltime_accounting_t walltime_accounting)
1703 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1705 gmx_localtop_t *top;
1706 gmx_enerdata_t *enerd;
1708 gmx_global_stat_t gstat;
1710 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1711 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1712 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1713 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1714 real a, b, c, maxdelta, delta;
1715 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1716 real dgdx, dgdg, sq, yr, beta;
1721 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1723 int start, end, number_steps;
1725 int i, k, m, n, nfmax, gf, step;
1730 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1735 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).");
1738 n = 3*state_global->natoms;
1739 nmaxcorr = inputrec->nbfgscorr;
1741 /* Allocate memory */
1742 /* Use pointers to real so we dont have to loop over both atoms and
1743 * dimensions all the time...
1744 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1745 * that point to the same memory.
1758 snew(rho, nmaxcorr);
1759 snew(alpha, nmaxcorr);
1762 for (i = 0; i < nmaxcorr; i++)
1768 for (i = 0; i < nmaxcorr; i++)
1777 init_em(fplog, LBFGS, cr, inputrec,
1778 state_global, top_global, &ems, &top, &f,
1779 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
1780 vsite, constr, NULL,
1781 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1782 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1783 * so we free some memory again.
1788 xx = (real *)state_global->x;
1792 end = mdatoms->homenr;
1794 /* Print to log file */
1795 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1797 do_log = do_ene = do_x = do_f = TRUE;
1799 /* Max number of steps */
1800 number_steps = inputrec->nsteps;
1802 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1804 for (i = start; i < end; i++)
1806 if (mdatoms->cFREEZE)
1808 gf = mdatoms->cFREEZE[i];
1810 for (m = 0; m < DIM; m++)
1812 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1817 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1821 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1826 construct_vsites(vsite, state_global->x, 1, NULL,
1827 top->idef.iparams, top->idef.il,
1828 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1831 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1832 /* do_force always puts the charge groups in the box and shifts again
1833 * We do not unshift, so molecules are always whole
1836 ems.s.x = state_global->x;
1838 evaluate_energy(fplog, cr,
1839 top_global, &ems, top,
1840 inputrec, nrnb, wcycle, gstat,
1841 vsite, constr, fcd, graph, mdatoms, fr,
1842 mu_tot, enerd, vir, pres, -1, TRUE);
1847 /* Copy stuff to the energy bin for easy printing etc. */
1848 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1849 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
1850 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1852 print_ebin_header(fplog, step, step);
1853 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1854 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1858 /* This is the starting energy */
1859 Epot = enerd->term[F_EPOT];
1865 /* Set the initial step.
1866 * since it will be multiplied by the non-normalized search direction
1867 * vector (force vector the first time), we scale it by the
1868 * norm of the force.
1873 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1874 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1875 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1876 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrtNumAtoms);
1877 fprintf(stderr, "\n");
1878 /* and copy to the log file too... */
1879 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1880 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1881 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrtNumAtoms);
1882 fprintf(fplog, "\n");
1885 // Point is an index to the memory of search directions, where 0 is the first one.
1888 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1889 for (i = 0; i < n; i++)
1893 dx[point][i] = ff[i]; /* Initial search direction */
1901 // Stepsize will be modified during the search, and actually it is not critical
1902 // (the main efficiency in the algorithm comes from changing directions), but
1903 // we still need an initial value, so estimate it as the inverse of the norm
1904 // so we take small steps where the potential fluctuates a lot.
1905 stepsize = 1.0/fnorm;
1907 /* Start the loop over BFGS steps.
1908 * Each successful step is counted, and we continue until
1909 * we either converge or reach the max number of steps.
1914 /* Set the gradient from the force */
1916 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1919 /* Write coordinates if necessary */
1920 do_x = do_per_step(step, inputrec->nstxout);
1921 do_f = do_per_step(step, inputrec->nstfout);
1926 mdof_flags |= MDOF_X;
1931 mdof_flags |= MDOF_F;
1936 mdof_flags |= MDOF_IMD;
1939 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1940 top_global, step, (real)step, state_global, state_global, f);
1942 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1944 /* make s a pointer to current search direction - point=0 first time we get here */
1947 // calculate line gradient in position A
1948 for (gpa = 0, i = 0; i < n; i++)
1953 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1954 * relative change in coordinate is smaller than precision
1956 for (minstep = 0, i = 0; i < n; i++)
1966 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1968 if (stepsize < minstep)
1974 // Before taking any steps along the line, store the old position
1975 for (i = 0; i < n; i++)
1982 for (i = 0; i < n; i++)
1987 /* Take a step downhill.
1988 * In theory, we should find the actual minimum of the function in this
1989 * direction, somewhere along the line.
1990 * That is quite possible, but it turns out to take 5-10 function evaluations
1991 * for each line. However, we dont really need to find the exact minimum -
1992 * it is much better to start a new BFGS step in a modified direction as soon
1993 * as we are close to it. This will save a lot of energy evaluations.
1995 * In practice, we just try to take a single step.
1996 * If it worked (i.e. lowered the energy), we increase the stepsize but
1997 * continue straight to the next BFGS step without trying to find any minimum,
1998 * i.e. we change the search direction too. If the line was smooth, it is
1999 * likely we are in a smooth region, and then it makes sense to take longer
2000 * steps in the modified search direction too.
2002 * If it didn't work (higher energy), there must be a minimum somewhere between
2003 * the old position and the new one. Then we need to start by finding a lower
2004 * value before we change search direction. Since the energy was apparently
2005 * quite rough, we need to decrease the step size.
2007 * Due to the finite numerical accuracy, it turns out that it is a good idea
2008 * to accept a SMALL increase in energy, if the derivative is still downhill.
2009 * This leads to lower final energies in the tests I've done. / Erik
2012 // State "A" is the first position along the line.
2013 // reference position along line is initially zero
2017 // Check stepsize first. We do not allow displacements
2018 // larger than emstep.
2022 // Pick a new position C by adding stepsize to A.
2025 // Calculate what the largest change in any individual coordinate
2026 // would be (translation along line * gradient along line)
2028 for (i = 0; i < n; i++)
2031 if (delta > maxdelta)
2036 // If any displacement is larger than the stepsize limit, reduce the step
2037 if (maxdelta > inputrec->em_stepsize)
2042 while (maxdelta > inputrec->em_stepsize);
2044 // Take a trial step and move the coordinate array xc[] to position C
2045 for (i = 0; i < n; i++)
2047 xc[i] = lastx[i] + c*s[i];
2051 // Calculate energy for the trial step in position C
2052 ems.s.x = (rvec *)xc;
2054 evaluate_energy(fplog, cr,
2055 top_global, &ems, top,
2056 inputrec, nrnb, wcycle, gstat,
2057 vsite, constr, fcd, graph, mdatoms, fr,
2058 mu_tot, enerd, vir, pres, step, FALSE);
2061 // Calc line gradient in position C
2062 for (gpc = 0, i = 0; i < n; i++)
2064 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2066 /* Sum the gradient along the line across CPUs */
2069 gmx_sumd(1, &gpc, cr);
2072 // This is the max amount of increase in energy we tolerate.
2073 // By allowing VERY small changes (close to numerical precision) we
2074 // frequently find even better (lower) final energies.
2075 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
2077 // Accept the step if the energy is lower in the new position C (compared to A),
2078 // or if it is not significantly higher and the line derivative is still negative.
2079 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
2081 // Great, we found a better energy. We no longer try to alter the
2082 // stepsize, but simply accept this new better position. The we select a new
2083 // search direction instead, which will be much more efficient than continuing
2084 // to take smaller steps along a line. Set fnorm based on the new C position,
2085 // which will be used to update the stepsize to 1/fnorm further down.
2091 // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2092 // or higher than in point A. In this case it is pointless to move to point C,
2093 // so we will have to do more iterations along the same line to find a smaller
2094 // value in the interval [A=0.0,C].
2095 // Here, A is still 0.0, but that will change when we do a search in the interval
2096 // [0.0,C] below. That search we will do by interpolation or bisection rather
2097 // than with the stepsize, so no need to modify it. For the next search direction
2098 // it will be reset to 1/fnorm anyway.
2104 // OK, if we didn't find a lower value we will have to locate one now - there must
2105 // be one in the interval [a,c].
2106 // The same thing is valid here, though: Don't spend dozens of iterations to find
2107 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2108 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2109 // I also have a safeguard for potentially really pathological functions so we never
2110 // take more than 20 steps before we give up.
2111 // If we already found a lower value we just skip this step and continue to the update.
2115 // Select a new trial point B in the interval [A,C].
2116 // If the derivatives at points a & c have different sign we interpolate to zero,
2117 // otherwise just do a bisection since there might be multiple minima/maxima
2118 // inside the interval.
2119 if (gpa < 0 && gpc > 0)
2121 b = a + gpa*(a-c)/(gpc-gpa);
2128 /* safeguard if interpolation close to machine accuracy causes errors:
2129 * never go outside the interval
2131 if (b <= a || b >= c)
2136 // Take a trial step to point B
2137 for (i = 0; i < n; i++)
2139 xb[i] = lastx[i] + b*s[i];
2143 // Calculate energy for the trial step in point B
2144 ems.s.x = (rvec *)xb;
2146 evaluate_energy(fplog, cr,
2147 top_global, &ems, top,
2148 inputrec, nrnb, wcycle, gstat,
2149 vsite, constr, fcd, graph, mdatoms, fr,
2150 mu_tot, enerd, vir, pres, step, FALSE);
2154 // Calculate gradient in point B
2155 for (gpb = 0, i = 0; i < n; i++)
2157 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2160 /* Sum the gradient along the line across CPUs */
2163 gmx_sumd(1, &gpb, cr);
2166 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2167 // at the new point B, and rename the endpoints of this new interval A and C.
2170 /* Replace c endpoint with b */
2174 /* swap coord pointers b/c */
2184 /* Replace a endpoint with b */
2188 /* swap coord pointers a/b */
2198 * Stop search as soon as we find a value smaller than the endpoints,
2199 * or if the tolerance is below machine precision.
2200 * Never run more than 20 steps, no matter what.
2204 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2206 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2208 /* OK. We couldn't find a significantly lower energy.
2209 * If ncorr==0 this was steepest descent, and then we give up.
2210 * If not, reset memory to restart as steepest descent before quitting.
2222 /* Search in gradient direction */
2223 for (i = 0; i < n; i++)
2225 dx[point][i] = ff[i];
2227 /* Reset stepsize */
2228 stepsize = 1.0/fnorm;
2233 /* Select min energy state of A & C, put the best in xx/ff/Epot
2239 for (i = 0; i < n; i++)
2250 for (i = 0; i < n; i++)
2264 for (i = 0; i < n; i++)
2272 /* Update the memory information, and calculate a new
2273 * approximation of the inverse hessian
2276 /* Have new data in Epot, xx, ff */
2277 if (ncorr < nmaxcorr)
2282 for (i = 0; i < n; i++)
2284 dg[point][i] = lastf[i]-ff[i];
2285 dx[point][i] *= step_taken;
2290 for (i = 0; i < n; i++)
2292 dgdg += dg[point][i]*dg[point][i];
2293 dgdx += dg[point][i]*dx[point][i];
2298 rho[point] = 1.0/dgdx;
2301 if (point >= nmaxcorr)
2307 for (i = 0; i < n; i++)
2314 /* Recursive update. First go back over the memory points */
2315 for (k = 0; k < ncorr; k++)
2324 for (i = 0; i < n; i++)
2326 sq += dx[cp][i]*p[i];
2329 alpha[cp] = rho[cp]*sq;
2331 for (i = 0; i < n; i++)
2333 p[i] -= alpha[cp]*dg[cp][i];
2337 for (i = 0; i < n; i++)
2342 /* And then go forward again */
2343 for (k = 0; k < ncorr; k++)
2346 for (i = 0; i < n; i++)
2348 yr += p[i]*dg[cp][i];
2352 beta = alpha[cp]-beta;
2354 for (i = 0; i < n; i++)
2356 p[i] += beta*dx[cp][i];
2366 for (i = 0; i < n; i++)
2370 dx[point][i] = p[i];
2378 /* Test whether the convergence criterion is met */
2379 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2381 /* Print it if necessary */
2386 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2387 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2388 step, Epot, fnorm/sqrtNumAtoms, fmax, nfmax+1);
2391 /* Store the new (lower) energies */
2392 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2393 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
2394 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2395 do_log = do_per_step(step, inputrec->nstlog);
2396 do_ene = do_per_step(step, inputrec->nstenergy);
2399 print_ebin_header(fplog, step, step);
2401 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2402 do_log ? fplog : NULL, step, step, eprNORMAL,
2403 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2406 /* Send x and E to IMD client, if bIMD is TRUE. */
2407 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2409 IMD_send_positions(inputrec->imd);
2412 // Reset stepsize in we are doing more iterations
2413 stepsize = 1.0/fnorm;
2415 /* Stop when the maximum force lies below tolerance.
2416 * If we have reached machine precision, converged is already set to true.
2418 converged = converged || (fmax < inputrec->em_tol);
2420 } /* End of the loop */
2422 /* IMD cleanup, if bIMD is TRUE. */
2423 IMD_finalize(inputrec->bIMD, inputrec->imd);
2427 step--; /* we never took that last step in this case */
2430 if (fmax > inputrec->em_tol)
2434 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2435 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2440 /* If we printed energy and/or logfile last step (which was the last step)
2441 * we don't have to do it again, but otherwise print the final values.
2443 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2445 print_ebin_header(fplog, step, step);
2447 if (!do_ene || !do_log) /* Write final energy file entries */
2449 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2450 !do_log ? fplog : NULL, step, step, eprNORMAL,
2451 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2454 /* Print some stuff... */
2457 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2461 * For accurate normal mode calculation it is imperative that we
2462 * store the last conformation into the full precision binary trajectory.
2464 * However, we should only do it if we did NOT already write this step
2465 * above (which we did if do_x or do_f was true).
2467 do_x = !do_per_step(step, inputrec->nstxout);
2468 do_f = !do_per_step(step, inputrec->nstfout);
2469 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2470 top_global, inputrec, step,
2471 &ems, state_global);
2475 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2476 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2477 number_steps, Epot, fmax, nfmax, fnorm/sqrtNumAtoms);
2478 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2479 number_steps, Epot, fmax, nfmax, fnorm/sqrtNumAtoms);
2481 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2484 finish_em(cr, outf, walltime_accounting, wcycle);
2486 /* To print the actual number of steps we needed somewhere */
2487 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2490 } /* That's all folks */
2492 /*! \brief Do steepest descents minimization
2493 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2494 int nfile, const t_filenm fnm[],
2495 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2497 gmx_vsite_t *vsite, gmx_constr_t constr,
2499 t_inputrec *inputrec,
2500 gmx_mtop_t *top_global, t_fcdata *fcd,
2501 t_state *state_global,
2503 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2506 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2507 real cpt_period, real max_hours,
2509 unsigned long Flags,
2510 gmx_walltime_accounting_t walltime_accounting)
2512 double do_steep(FILE *fplog, t_commrec *cr, const gmx::MDLogger gmx_unused &mdlog,
2513 int nfile, const t_filenm fnm[],
2514 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2515 int gmx_unused nstglobalcomm,
2516 gmx_vsite_t *vsite, gmx_constr_t constr,
2517 int gmx_unused stepout,
2518 t_inputrec *inputrec,
2519 gmx_mtop_t *top_global, t_fcdata *fcd,
2520 t_state *state_global,
2522 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2523 gmx_edsam_t gmx_unused ed,
2525 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2526 gmx_membed_t gmx_unused *membed,
2527 real gmx_unused cpt_period, real gmx_unused max_hours,
2529 unsigned long gmx_unused Flags,
2530 gmx_walltime_accounting_t walltime_accounting)
2532 const char *SD = "Steepest Descents";
2533 em_state_t *s_min, *s_try;
2534 gmx_localtop_t *top;
2535 gmx_enerdata_t *enerd;
2537 gmx_global_stat_t gstat;
2543 gmx_bool bDone, bAbort, do_x, do_f;
2548 int steps_accepted = 0;
2550 s_min = init_em_state();
2551 s_try = init_em_state();
2553 /* Init em and store the local state in s_try */
2554 init_em(fplog, SD, cr, inputrec,
2555 state_global, top_global, s_try, &top, &f,
2556 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
2557 vsite, constr, NULL,
2558 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
2560 /* Print to log file */
2561 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2563 /* Set variables for stepsize (in nm). This is the largest
2564 * step that we are going to make in any direction.
2566 ustep = inputrec->em_stepsize;
2569 /* Max number of steps */
2570 nsteps = inputrec->nsteps;
2574 /* Print to the screen */
2575 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2579 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2582 /**** HERE STARTS THE LOOP ****
2583 * count is the counter for the number of steps
2584 * bDone will be TRUE when the minimization has converged
2585 * bAbort will be TRUE when nsteps steps have been performed or when
2586 * the stepsize becomes smaller than is reasonable for machine precision
2591 while (!bDone && !bAbort)
2593 bAbort = (nsteps >= 0) && (count == nsteps);
2595 /* set new coordinates, except for first step */
2596 bool validStep = true;
2600 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2601 s_min, stepsize, s_min->f, s_try,
2602 constr, top, nrnb, wcycle, count);
2607 evaluate_energy(fplog, cr,
2608 top_global, s_try, top,
2609 inputrec, nrnb, wcycle, gstat,
2610 vsite, constr, fcd, graph, mdatoms, fr,
2611 mu_tot, enerd, vir, pres, count, count == 0);
2615 // Signal constraint error during stepping with energy=inf
2616 s_try->epot = std::numeric_limits<real>::infinity();
2621 print_ebin_header(fplog, count, count);
2626 s_min->epot = s_try->epot;
2629 /* Print it if necessary */
2634 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2635 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2636 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2640 if ( (count == 0) || (s_try->epot < s_min->epot) )
2642 /* Store the new (lower) energies */
2643 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2644 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2645 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2647 /* Prepare IMD energy record, if bIMD is TRUE. */
2648 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2650 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2651 do_per_step(steps_accepted, inputrec->nstdisreout),
2652 do_per_step(steps_accepted, inputrec->nstorireout),
2653 fplog, count, count, eprNORMAL,
2654 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2659 /* Now if the new energy is smaller than the previous...
2660 * or if this is the first step!
2661 * or if we did random steps!
2664 if ( (count == 0) || (s_try->epot < s_min->epot) )
2668 /* Test whether the convergence criterion is met... */
2669 bDone = (s_try->fmax < inputrec->em_tol);
2671 /* Copy the arrays for force, positions and energy */
2672 /* The 'Min' array always holds the coords and forces of the minimal
2674 swap_em_state(s_min, s_try);
2680 /* Write to trn, if necessary */
2681 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2682 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2683 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2684 top_global, inputrec, count,
2685 s_min, state_global);
2689 /* If energy is not smaller make the step smaller... */
2692 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2694 /* Reload the old state */
2695 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2696 s_min, top, mdatoms, fr, vsite, constr,
2701 /* Determine new step */
2702 stepsize = ustep/s_min->fmax;
2704 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2706 if (count == nsteps || ustep < 1e-12)
2708 if (count == nsteps || ustep < 1e-6)
2713 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2714 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2719 /* Send IMD energies and positions, if bIMD is TRUE. */
2720 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2722 IMD_send_positions(inputrec->imd);
2726 } /* End of the loop */
2728 /* IMD cleanup, if bIMD is TRUE. */
2729 IMD_finalize(inputrec->bIMD, inputrec->imd);
2731 /* Print some data... */
2734 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2736 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2737 top_global, inputrec, count,
2738 s_min, state_global);
2742 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2743 fnormn = s_min->fnorm/sqrtNumAtoms;
2745 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2746 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2747 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2748 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2751 finish_em(cr, outf, walltime_accounting, wcycle);
2753 /* To print the actual number of steps we needed somewhere */
2754 inputrec->nsteps = count;
2756 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2759 } /* That's all folks */
2761 /*! \brief Do normal modes analysis
2762 \copydoc integrator_t(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2763 int nfile, const t_filenm fnm[],
2764 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2766 gmx_vsite_t *vsite, gmx_constr_t constr,
2768 t_inputrec *inputrec,
2769 gmx_mtop_t *top_global, t_fcdata *fcd,
2770 t_state *state_global,
2772 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2775 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2776 real cpt_period, real max_hours,
2778 unsigned long Flags,
2779 gmx_walltime_accounting_t walltime_accounting)
2781 double do_nm(FILE *fplog, t_commrec *cr, const gmx::MDLogger &mdlog,
2782 int nfile, const t_filenm fnm[],
2783 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2784 int gmx_unused nstglobalcomm,
2785 gmx_vsite_t *vsite, gmx_constr_t constr,
2786 int gmx_unused stepout,
2787 t_inputrec *inputrec,
2788 gmx_mtop_t *top_global, t_fcdata *fcd,
2789 t_state *state_global,
2791 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2792 gmx_edsam_t gmx_unused ed,
2794 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2795 gmx_membed_t gmx_unused *membed,
2796 real gmx_unused cpt_period, real gmx_unused max_hours,
2798 unsigned long gmx_unused Flags,
2799 gmx_walltime_accounting_t walltime_accounting)
2801 const char *NM = "Normal Mode Analysis";
2804 gmx_localtop_t *top;
2805 gmx_enerdata_t *enerd;
2807 gmx_global_stat_t gstat;
2812 gmx_bool bSparse; /* use sparse matrix storage format */
2814 gmx_sparsematrix_t * sparse_matrix = NULL;
2815 real * full_matrix = NULL;
2816 em_state_t * state_work;
2818 /* added with respect to mdrun */
2820 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2822 bool bIsMaster = MASTER(cr);
2826 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2829 state_work = init_em_state();
2831 gmx_shellfc_t *shellfc;
2833 /* Init em and store the local state in state_minimum */
2834 init_em(fplog, NM, cr, inputrec,
2835 state_global, top_global, state_work, &top,
2837 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat,
2838 vsite, constr, &shellfc,
2839 nfile, fnm, &outf, NULL, imdport, Flags, wcycle);
2841 std::vector<size_t> atom_index = get_atom_index(top_global);
2842 snew(fneg, atom_index.size());
2843 snew(dfdx, atom_index.size());
2849 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2850 " which MIGHT not be accurate enough for normal mode analysis.\n"
2851 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2852 " are fairly modest even if you recompile in double precision.\n\n");
2856 /* Check if we can/should use sparse storage format.
2858 * Sparse format is only useful when the Hessian itself is sparse, which it
2859 * will be when we use a cutoff.
2860 * For small systems (n<1000) it is easier to always use full matrix format, though.
2862 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2864 GMX_LOG(mdlog.warning).appendText("Non-cutoff electrostatics used, forcing full Hessian format.");
2867 else if (atom_index.size() < 1000)
2869 GMX_LOG(mdlog.warning).appendTextFormatted("Small system size (N=%d), using full Hessian format.",
2875 GMX_LOG(mdlog.warning).appendText("Using compressed symmetric sparse Hessian format.");
2879 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2880 sz = DIM*atom_index.size();
2882 fprintf(stderr, "Allocating Hessian memory...\n\n");
2886 sparse_matrix = gmx_sparsematrix_init(sz);
2887 sparse_matrix->compressed_symmetric = TRUE;
2891 snew(full_matrix, sz*sz);
2898 /* Write start time and temperature */
2899 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2901 /* fudge nr of steps to nr of atoms */
2902 inputrec->nsteps = atom_index.size()*2;
2906 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2907 *(top_global->name), (int)inputrec->nsteps);
2910 nnodes = cr->nnodes;
2912 /* Make evaluate_energy do a single node force calculation */
2914 evaluate_energy(fplog, cr,
2915 top_global, state_work, top,
2916 inputrec, nrnb, wcycle, gstat,
2917 vsite, constr, fcd, graph, mdatoms, fr,
2918 mu_tot, enerd, vir, pres, -1, TRUE);
2919 cr->nnodes = nnodes;
2921 /* if forces are not small, warn user */
2922 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2924 GMX_LOG(mdlog.warning).appendTextFormatted("Maximum force:%12.5e", state_work->fmax);
2925 if (state_work->fmax > 1.0e-3)
2927 GMX_LOG(mdlog.warning).appendText(
2928 "The force is probably not small enough to "
2929 "ensure that you are at a minimum.\n"
2930 "Be aware that negative eigenvalues may occur\n"
2931 "when the resulting matrix is diagonalized.");
2934 /***********************************************************
2936 * Loop over all pairs in matrix
2938 * do_force called twice. Once with positive and
2939 * once with negative displacement
2941 ************************************************************/
2943 /* Steps are divided one by one over the nodes */
2945 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2947 size_t atom = atom_index[aid];
2948 for (size_t d = 0; d < DIM; d++)
2950 gmx_bool bBornRadii = FALSE;
2951 gmx_int64_t step = 0;
2952 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2955 x_min = state_work->s.x[atom][d];
2957 for (unsigned int dx = 0; (dx < 2); dx++)
2961 state_work->s.x[atom][d] = x_min - der_range;
2965 state_work->s.x[atom][d] = x_min + der_range;
2968 /* Make evaluate_energy do a single node force calculation */
2972 /* Now is the time to relax the shells */
2973 (void) relax_shell_flexcon(fplog, cr, bVerbose, step,
2974 inputrec, bNS, force_flags,
2977 &state_work->s, state_work->f, vir, mdatoms,
2978 nrnb, wcycle, graph, &top_global->groups,
2979 shellfc, fr, bBornRadii, t, mu_tot,
2986 evaluate_energy(fplog, cr,
2987 top_global, state_work, top,
2988 inputrec, nrnb, wcycle, gstat,
2989 vsite, constr, fcd, graph, mdatoms, fr,
2990 mu_tot, enerd, vir, pres, atom*2+dx, FALSE);
2993 cr->nnodes = nnodes;
2997 for (size_t i = 0; i < atom_index.size(); i++)
2999 copy_rvec(state_work->f[atom_index[i]], fneg[i]);
3004 /* x is restored to original */
3005 state_work->s.x[atom][d] = x_min;
3007 for (size_t j = 0; j < atom_index.size(); j++)
3009 for (size_t k = 0; (k < DIM); k++)
3012 -(state_work->f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
3019 #define mpi_type GMX_MPI_REAL
3020 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
3021 cr->nodeid, cr->mpi_comm_mygroup);
3026 for (node = 0; (node < nnodes && atom+node < atom_index.size()); node++)
3032 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
3033 cr->mpi_comm_mygroup, &stat);
3038 row = (atom + node)*DIM + d;
3040 for (size_t j = 0; j < atom_index.size(); j++)
3042 for (size_t k = 0; k < DIM; k++)
3048 if (col >= row && dfdx[j][k] != 0.0)
3050 gmx_sparsematrix_increment_value(sparse_matrix,
3051 row, col, dfdx[j][k]);
3056 full_matrix[row*sz+col] = dfdx[j][k];
3063 if (bVerbose && fplog)
3068 /* write progress */
3069 if (bIsMaster && bVerbose)
3071 fprintf(stderr, "\rFinished step %d out of %d",
3072 static_cast<int>(std::min(atom+nnodes, atom_index.size())),
3073 static_cast<int>(atom_index.size()));
3080 fprintf(stderr, "\n\nWriting Hessian...\n");
3081 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
3084 finish_em(cr, outf, walltime_accounting, wcycle);
3086 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);