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/md_logging.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gmxlib/nrnb.h"
67 #include "gromacs/imd/imd.h"
68 #include "gromacs/linearalgebra/sparsematrix.h"
69 #include "gromacs/listed-forces/manage-threading.h"
70 #include "gromacs/math/functions.h"
71 #include "gromacs/math/vec.h"
72 #include "gromacs/mdlib/constr.h"
73 #include "gromacs/mdlib/force.h"
74 #include "gromacs/mdlib/forcerec.h"
75 #include "gromacs/mdlib/gmx_omp_nthreads.h"
76 #include "gromacs/mdlib/md_support.h"
77 #include "gromacs/mdlib/mdatoms.h"
78 #include "gromacs/mdlib/mdebin.h"
79 #include "gromacs/mdlib/mdrun.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/utility/cstringutil.h"
96 #include "gromacs/utility/exceptions.h"
97 #include "gromacs/utility/fatalerror.h"
98 #include "gromacs/utility/smalloc.h"
100 //! Utility structure for manipulating states during EM
102 //! Copy of the global state
108 //! Norm of the force
116 //! Initiate em_state_t structure and return pointer to it
117 static em_state_t *init_em_state()
123 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
124 snew(ems->s.lambda, efptNR);
129 //! Print the EM starting conditions
130 static void print_em_start(FILE *fplog,
132 gmx_walltime_accounting_t walltime_accounting,
133 gmx_wallcycle_t wcycle,
136 walltime_accounting_start(walltime_accounting);
137 wallcycle_start(wcycle, ewcRUN);
138 print_start(fplog, cr, walltime_accounting, name);
141 //! Stop counting time for EM
142 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
143 gmx_wallcycle_t wcycle)
145 wallcycle_stop(wcycle, ewcRUN);
147 walltime_accounting_end(walltime_accounting);
150 //! Printing a log file and console header
151 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
154 fprintf(out, "%s:\n", minimizer);
155 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
156 fprintf(out, " Number of steps = %12d\n", nsteps);
159 //! Print warning message
160 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
166 "\nEnergy minimization reached the maximum number "
167 "of steps before the forces reached the requested "
168 "precision Fmax < %g.\n", ftol);
173 "\nEnergy minimization has stopped, but the forces have "
174 "not converged to the requested precision Fmax < %g (which "
175 "may not be possible for your system). It stopped "
176 "because the algorithm tried to make a new step whose size "
177 "was too small, or there was no change in the energy since "
178 "last step. Either way, we regard the minimization as "
179 "converged to within the available machine precision, "
180 "given your starting configuration and EM parameters.\n%s%s",
182 sizeof(real) < sizeof(double) ?
183 "\nDouble precision normally gives you higher accuracy, but "
184 "this is often not needed for preparing to run molecular "
188 "You might need to increase your constraint accuracy, or turn\n"
189 "off constraints altogether (set constraints = none in mdp file)\n" :
192 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
195 //! Print message about convergence of the EM
196 static void print_converged(FILE *fp, const char *alg, real ftol,
197 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
198 real epot, real fmax, int nfmax, real fnorm)
200 char buf[STEPSTRSIZE];
204 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
205 alg, ftol, gmx_step_str(count, buf));
207 else if (count < nsteps)
209 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
210 "but did not reach the requested Fmax < %g.\n",
211 alg, gmx_step_str(count, buf), ftol);
215 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
216 alg, ftol, gmx_step_str(count, buf));
220 fprintf(fp, "Potential Energy = %21.14e\n", epot);
221 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
222 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
224 fprintf(fp, "Potential Energy = %14.7e\n", epot);
225 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
226 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
230 //! Compute the norm and max of the force array in parallel
231 static void get_f_norm_max(t_commrec *cr,
232 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
233 real *fnorm, real *fmax, int *a_fmax)
237 int la_max, a_max, start, end, i, m, gf;
239 /* This routine finds the largest force and returns it.
240 * On parallel machines the global max is taken.
246 end = mdatoms->homenr;
247 if (mdatoms->cFREEZE)
249 for (i = start; i < end; i++)
251 gf = mdatoms->cFREEZE[i];
253 for (m = 0; m < DIM; m++)
255 if (!opts->nFreeze[gf][m])
257 fam += gmx::square(f[i][m]);
270 for (i = start; i < end; i++)
282 if (la_max >= 0 && DOMAINDECOMP(cr))
284 a_max = cr->dd->gatindex[la_max];
292 snew(sum, 2*cr->nnodes+1);
293 sum[2*cr->nodeid] = fmax2;
294 sum[2*cr->nodeid+1] = a_max;
295 sum[2*cr->nnodes] = fnorm2;
296 gmx_sumd(2*cr->nnodes+1, sum, cr);
297 fnorm2 = sum[2*cr->nnodes];
298 /* Determine the global maximum */
299 for (i = 0; i < cr->nnodes; i++)
301 if (sum[2*i] > fmax2)
304 a_max = (int)(sum[2*i+1] + 0.5);
312 *fnorm = sqrt(fnorm2);
324 //! Compute the norm of the force
325 static void get_state_f_norm_max(t_commrec *cr,
326 t_grpopts *opts, t_mdatoms *mdatoms,
329 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
332 //! Initialize the energy minimization
333 void init_em(FILE *fplog, const char *title,
334 t_commrec *cr, t_inputrec *ir,
335 t_state *state_global, gmx_mtop_t *top_global,
336 em_state_t *ems, gmx_localtop_t **top,
338 t_nrnb *nrnb, rvec mu_tot,
339 t_forcerec *fr, gmx_enerdata_t **enerd,
340 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
341 gmx_vsite_t *vsite, gmx_constr_t constr,
342 int nfile, const t_filenm fnm[],
343 gmx_mdoutf_t *outf, t_mdebin **mdebin,
344 int imdport, unsigned long gmx_unused Flags,
345 gmx_wallcycle_t wcycle)
352 fprintf(fplog, "Initiating %s\n", title);
355 state_global->ngtc = 0;
357 /* Initialize lambda variables */
358 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
362 /* Interactive molecular dynamics */
363 init_IMD(ir, cr, top_global, fplog, 1, state_global->x,
364 nfile, fnm, NULL, imdport, Flags);
366 if (DOMAINDECOMP(cr))
368 *top = dd_init_local_top(top_global);
370 dd_init_local_state(cr->dd, state_global, &ems->s);
374 /* Distribute the charge groups over the nodes from the master node */
375 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
376 state_global, top_global, ir,
377 &ems->s, &ems->f, mdatoms, *top,
380 dd_store_state(cr->dd, &ems->s);
386 snew(*f, top_global->natoms);
388 /* Just copy the state */
389 ems->s = *state_global;
390 /* We need to allocate one element extra, since we might use
391 * (unaligned) 4-wide SIMD loads to access rvec entries.
393 snew(ems->s.x, ems->s.nalloc + 1);
394 snew(ems->f, ems->s.nalloc+1);
395 snew(ems->s.v, ems->s.nalloc+1);
396 for (i = 0; i < state_global->natoms; i++)
398 copy_rvec(state_global->x[i], ems->s.x[i]);
400 copy_mat(state_global->box, ems->s.box);
402 *top = gmx_mtop_generate_local_top(top_global, ir->efep != efepNO);
404 setup_bonded_threading(fr, &(*top)->idef);
406 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
408 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
415 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
416 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
420 set_vsite_top(vsite, *top, mdatoms, cr);
426 if (ir->eConstrAlg == econtSHAKE &&
427 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
429 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
430 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
433 if (!DOMAINDECOMP(cr))
435 set_constraints(constr, *top, ir, mdatoms, cr);
438 if (!ir->bContinuation)
440 /* Constrain the starting coordinates */
442 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
443 ir, cr, -1, 0, 1.0, mdatoms,
444 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
445 ems->s.lambda[efptFEP], &dvdl_constr,
446 NULL, NULL, nrnb, econqCoord);
452 *gstat = global_stat_init(ir);
459 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, NULL, wcycle);
462 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
467 /* Init bin for energy stuff */
468 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
472 calc_shifts(ems->s.box, fr->shift_vec);
475 //! Finalize the minimization
476 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
477 gmx_walltime_accounting_t walltime_accounting,
478 gmx_wallcycle_t wcycle)
480 if (!(cr->duty & DUTY_PME))
482 /* Tell the PME only node to finish */
483 gmx_pme_send_finish(cr);
488 em_time_end(walltime_accounting, wcycle);
491 //! Swap two different EM states during minimization
492 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
501 //! Copy coordinate from an EM state to a "normal" state structure
502 static void copy_em_coords(em_state_t *ems, t_state *state)
506 for (i = 0; (i < state->natoms); i++)
508 copy_rvec(ems->s.x[i], state->x[i]);
512 //! Save the EM trajectory
513 static void write_em_traj(FILE *fplog, t_commrec *cr,
515 gmx_bool bX, gmx_bool bF, const char *confout,
516 gmx_mtop_t *top_global,
517 t_inputrec *ir, gmx_int64_t step,
519 t_state *state_global)
522 gmx_bool bIMDout = FALSE;
525 /* Shall we do IMD output? */
528 bIMDout = do_per_step(step, IMD_get_step(ir->imd->setup));
531 if ((bX || bF || bIMDout || confout != NULL) && !DOMAINDECOMP(cr))
533 copy_em_coords(state, state_global);
539 mdof_flags |= MDOF_X;
543 mdof_flags |= MDOF_F;
546 /* If we want IMD output, set appropriate MDOF flag */
549 mdof_flags |= MDOF_IMD;
552 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
553 top_global, step, (double)step,
554 &state->s, state_global, state->f);
556 if (confout != NULL && MASTER(cr))
558 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
560 /* Make molecules whole only for confout writing */
561 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
565 write_sto_conf_mtop(confout,
566 *top_global->name, top_global,
567 state_global->x, NULL, ir->ePBC, state_global->box);
571 //! \brief Do one minimization step
573 // \returns true when the step succeeded, false when a constraint error occurred
574 static bool do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
576 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
577 gmx_constr_t constr, gmx_localtop_t *top,
578 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
585 int nthreads gmx_unused;
587 bool validStep = true;
592 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
594 gmx_incons("state mismatch in do_em_step");
597 s2->flags = s1->flags;
599 if (s2->nalloc != s1->nalloc)
601 s2->nalloc = s1->nalloc;
602 /* We need to allocate one element extra, since we might use
603 * (unaligned) 4-wide SIMD loads to access rvec entries.
605 srenew(s2->x, s1->nalloc + 1);
606 srenew(ems2->f, s1->nalloc);
607 if (s2->flags & (1<<estCGP))
609 srenew(s2->cg_p, s1->nalloc + 1);
613 s2->natoms = s1->natoms;
614 copy_mat(s1->box, s2->box);
615 /* Copy free energy state */
616 for (int i = 0; i < efptNR; i++)
618 s2->lambda[i] = s1->lambda[i];
620 copy_mat(s1->box, s2->box);
625 // cppcheck-suppress unreadVariable
626 nthreads = gmx_omp_nthreads_get(emntUpdate);
627 #pragma omp parallel num_threads(nthreads)
633 #pragma omp for schedule(static) nowait
634 for (int i = start; i < end; i++)
642 for (int m = 0; m < DIM; m++)
644 if (ir->opts.nFreeze[gf][m])
650 x2[i][m] = x1[i][m] + a*f[i][m];
654 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
657 if (s2->flags & (1<<estCGP))
659 /* Copy the CG p vector */
662 #pragma omp for schedule(static) nowait
663 for (int i = start; i < end; i++)
665 // Trivial OpenMP block that does not throw
666 copy_rvec(p1[i], p2[i]);
670 if (DOMAINDECOMP(cr))
672 s2->ddp_count = s1->ddp_count;
673 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
676 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
679 /* We need to allocate one element extra, since we might use
680 * (unaligned) 4-wide SIMD loads to access rvec entries.
682 srenew(s2->cg_gl, s2->cg_gl_nalloc + 1);
684 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
687 s2->ncg_gl = s1->ncg_gl;
688 #pragma omp for schedule(static) nowait
689 for (int i = 0; i < s2->ncg_gl; i++)
691 s2->cg_gl[i] = s1->cg_gl[i];
693 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
699 wallcycle_start(wcycle, ewcCONSTR);
702 constrain(NULL, TRUE, TRUE, constr, &top->idef,
703 ir, cr, count, 0, 1.0, md,
704 s1->x, s2->x, NULL, bMolPBC, s2->box,
705 s2->lambda[efptBONDED], &dvdl_constr,
706 NULL, NULL, nrnb, econqCoord);
707 wallcycle_stop(wcycle, ewcCONSTR);
709 // We should move this check to the different minimizers
710 if (!validStep && ir->eI != eiSteep)
712 gmx_fatal(FARGS, "The coordinates could not be constrained. Minimizer '%s' can not handle constraint failures, use minimizer '%s' before using '%s'.",
713 EI(ir->eI), EI(eiSteep), EI(ir->eI));
720 //! Prepare EM for using domain decomposition parallellization
721 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
722 gmx_mtop_t *top_global, t_inputrec *ir,
723 em_state_t *ems, gmx_localtop_t *top,
724 t_mdatoms *mdatoms, t_forcerec *fr,
725 gmx_vsite_t *vsite, gmx_constr_t constr,
726 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
728 /* Repartition the domain decomposition */
729 dd_partition_system(fplog, step, cr, FALSE, 1,
730 NULL, top_global, ir,
732 mdatoms, top, fr, vsite, constr,
733 nrnb, wcycle, FALSE);
734 dd_store_state(cr->dd, &ems->s);
737 //! De one energy evaluation
738 static void evaluate_energy(FILE *fplog, t_commrec *cr,
739 gmx_mtop_t *top_global,
740 em_state_t *ems, gmx_localtop_t *top,
741 t_inputrec *inputrec,
742 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
743 gmx_global_stat_t gstat,
744 gmx_vsite_t *vsite, gmx_constr_t constr,
746 t_graph *graph, t_mdatoms *mdatoms,
747 t_forcerec *fr, rvec mu_tot,
748 gmx_enerdata_t *enerd, tensor vir, tensor pres,
749 gmx_int64_t count, gmx_bool bFirst)
753 tensor force_vir, shake_vir, ekin;
754 real dvdl_constr, prescorr, enercorr, dvdlcorr;
757 /* Set the time to the initial time, the time does not change during EM */
758 t = inputrec->init_t;
761 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
763 /* This is the first state or an old state used before the last ns */
769 if (inputrec->nstlist > 0)
777 construct_vsites(vsite, ems->s.x, 1, NULL,
778 top->idef.iparams, top->idef.il,
779 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
782 if (DOMAINDECOMP(cr) && bNS)
784 /* Repartition the domain decomposition */
785 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
786 ems, top, mdatoms, fr, vsite, constr,
790 /* Calc force & energy on new trial position */
791 /* do_force always puts the charge groups in the box and shifts again
792 * We do not unshift, so molecules are always whole in congrad.c
794 do_force(fplog, cr, inputrec,
795 count, nrnb, wcycle, top, &top_global->groups,
796 ems->s.box, ems->s.x, &ems->s.hist,
797 ems->f, force_vir, mdatoms, enerd, fcd,
798 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
799 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
800 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
801 (bNS ? GMX_FORCE_NS : 0));
803 /* Clear the unused shake virial and pressure */
804 clear_mat(shake_vir);
807 /* Communicate stuff when parallel */
808 if (PAR(cr) && inputrec->eI != eiNM)
810 wallcycle_start(wcycle, ewcMoveE);
812 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
813 inputrec, NULL, NULL, NULL, 1, &terminate,
819 wallcycle_stop(wcycle, ewcMoveE);
822 /* Calculate long range corrections to pressure and energy */
823 calc_dispcorr(inputrec, fr, ems->s.box, ems->s.lambda[efptVDW],
824 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
825 enerd->term[F_DISPCORR] = enercorr;
826 enerd->term[F_EPOT] += enercorr;
827 enerd->term[F_PRES] += prescorr;
828 enerd->term[F_DVDL] += dvdlcorr;
830 ems->epot = enerd->term[F_EPOT];
834 /* Project out the constraint components of the force */
835 wallcycle_start(wcycle, ewcCONSTR);
837 constrain(NULL, FALSE, FALSE, constr, &top->idef,
838 inputrec, cr, count, 0, 1.0, mdatoms,
839 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
840 ems->s.lambda[efptBONDED], &dvdl_constr,
841 NULL, &shake_vir, nrnb, econqForceDispl);
842 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
843 m_add(force_vir, shake_vir, vir);
844 wallcycle_stop(wcycle, ewcCONSTR);
848 copy_mat(force_vir, vir);
852 enerd->term[F_PRES] =
853 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
855 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
857 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
859 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
863 //! Parallel utility summing energies and forces
864 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
865 gmx_mtop_t *top_global,
866 em_state_t *s_min, em_state_t *s_b)
870 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
872 unsigned char *grpnrFREEZE;
876 fprintf(debug, "Doing reorder_partsum\n");
882 cgs_gl = dd_charge_groups_global(cr->dd);
883 index = cgs_gl->index;
885 /* Collect fm in a global vector fmg.
886 * This conflicts with the spirit of domain decomposition,
887 * but to fully optimize this a much more complicated algorithm is required.
889 snew(fmg, top_global->natoms);
891 ncg = s_min->s.ncg_gl;
892 cg_gl = s_min->s.cg_gl;
894 for (c = 0; c < ncg; c++)
899 for (a = a0; a < a1; a++)
901 copy_rvec(fm[i], fmg[a]);
905 gmx_sum(top_global->natoms*3, fmg[0], cr);
907 /* Now we will determine the part of the sum for the cgs in state s_b */
909 cg_gl = s_b->s.cg_gl;
913 grpnrFREEZE = top_global->groups.grpnr[egcFREEZE];
914 for (c = 0; c < ncg; c++)
919 for (a = a0; a < a1; a++)
921 if (mdatoms->cFREEZE && grpnrFREEZE)
925 for (m = 0; m < DIM; m++)
927 if (!opts->nFreeze[gf][m])
929 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
941 //! Print some stuff, like beta, whatever that means.
942 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
943 gmx_mtop_t *top_global,
944 em_state_t *s_min, em_state_t *s_b)
950 /* This is just the classical Polak-Ribiere calculation of beta;
951 * it looks a bit complicated since we take freeze groups into account,
952 * and might have to sum it in parallel runs.
955 if (!DOMAINDECOMP(cr) ||
956 (s_min->s.ddp_count == cr->dd->ddp_count &&
957 s_b->s.ddp_count == cr->dd->ddp_count))
963 /* This part of code can be incorrect with DD,
964 * since the atom ordering in s_b and s_min might differ.
966 for (i = 0; i < mdatoms->homenr; i++)
968 if (mdatoms->cFREEZE)
970 gf = mdatoms->cFREEZE[i];
972 for (m = 0; m < DIM; m++)
974 if (!opts->nFreeze[gf][m])
976 sum += (fb[i][m] - fm[i][m])*fb[i][m];
983 /* We need to reorder cgs while summing */
984 sum = reorder_partsum(cr, opts, mdatoms, top_global, s_min, s_b);
988 gmx_sumd(1, &sum, cr);
991 return sum/gmx::square(s_min->fnorm);
997 /*! \brief Do conjugate gradients minimization
998 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
999 int nfile, const t_filenm fnm[],
1000 const gmx_output_env_t *oenv, gmx_bool bVerbose,
1002 gmx_vsite_t *vsite, gmx_constr_t constr,
1004 t_inputrec *inputrec,
1005 gmx_mtop_t *top_global, t_fcdata *fcd,
1006 t_state *state_global,
1008 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1011 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1012 gmx_membed_t gmx_unused *membed,
1013 real cpt_period, real max_hours,
1015 unsigned long Flags,
1016 gmx_walltime_accounting_t walltime_accounting)
1018 double do_cg(FILE *fplog, t_commrec *cr,
1019 int nfile, const t_filenm fnm[],
1020 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
1021 int gmx_unused nstglobalcomm,
1022 gmx_vsite_t *vsite, gmx_constr_t constr,
1023 int gmx_unused stepout,
1024 t_inputrec *inputrec,
1025 gmx_mtop_t *top_global, t_fcdata *fcd,
1026 t_state *state_global,
1028 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1029 gmx_edsam_t gmx_unused ed,
1031 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1032 gmx_membed_t gmx_unused *membed,
1033 real gmx_unused cpt_period, real gmx_unused max_hours,
1035 unsigned long gmx_unused Flags,
1036 gmx_walltime_accounting_t walltime_accounting)
1038 const char *CG = "Polak-Ribiere Conjugate Gradients";
1040 em_state_t *s_min, *s_a, *s_b, *s_c;
1041 gmx_localtop_t *top;
1042 gmx_enerdata_t *enerd;
1044 gmx_global_stat_t gstat;
1047 double gpa, gpb, gpc, tmp, minstep;
1050 real a, b, c, beta = 0.0;
1054 gmx_bool converged, foundlower;
1056 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
1058 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
1060 int i, m, gf, step, nminstep;
1064 s_min = init_em_state();
1065 s_a = init_em_state();
1066 s_b = init_em_state();
1067 s_c = init_em_state();
1069 /* Init em and store the local state in s_min */
1070 init_em(fplog, CG, cr, inputrec,
1071 state_global, top_global, s_min, &top, &f,
1072 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1073 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1075 /* Print to log file */
1076 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1078 /* Max number of steps */
1079 number_steps = inputrec->nsteps;
1083 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1087 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1090 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1091 /* do_force always puts the charge groups in the box and shifts again
1092 * We do not unshift, so molecules are always whole in congrad.c
1094 evaluate_energy(fplog, cr,
1095 top_global, s_min, top,
1096 inputrec, nrnb, wcycle, gstat,
1097 vsite, constr, fcd, graph, mdatoms, fr,
1098 mu_tot, enerd, vir, pres, -1, TRUE);
1103 /* Copy stuff to the energy bin for easy printing etc. */
1104 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1105 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1106 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1108 print_ebin_header(fplog, step, step);
1109 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1110 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1114 /* Estimate/guess the initial stepsize */
1115 stepsize = inputrec->em_stepsize/s_min->fnorm;
1119 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1120 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1121 s_min->fmax, s_min->a_fmax+1);
1122 fprintf(stderr, " F-Norm = %12.5e\n",
1123 s_min->fnorm/sqrtNumAtoms);
1124 fprintf(stderr, "\n");
1125 /* and copy to the log file too... */
1126 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1127 s_min->fmax, s_min->a_fmax+1);
1128 fprintf(fplog, " F-Norm = %12.5e\n",
1129 s_min->fnorm/sqrtNumAtoms);
1130 fprintf(fplog, "\n");
1132 /* Start the loop over CG steps.
1133 * Each successful step is counted, and we continue until
1134 * we either converge or reach the max number of steps.
1137 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1140 /* start taking steps in a new direction
1141 * First time we enter the routine, beta=0, and the direction is
1142 * simply the negative gradient.
1145 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1150 for (i = 0; i < mdatoms->homenr; i++)
1152 if (mdatoms->cFREEZE)
1154 gf = mdatoms->cFREEZE[i];
1156 for (m = 0; m < DIM; m++)
1158 if (!inputrec->opts.nFreeze[gf][m])
1160 p[i][m] = sf[i][m] + beta*p[i][m];
1161 gpa -= p[i][m]*sf[i][m];
1162 /* f is negative gradient, thus the sign */
1171 /* Sum the gradient along the line across CPUs */
1174 gmx_sumd(1, &gpa, cr);
1177 /* Calculate the norm of the search vector */
1178 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1180 /* Just in case stepsize reaches zero due to numerical precision... */
1183 stepsize = inputrec->em_stepsize/pnorm;
1187 * Double check the value of the derivative in the search direction.
1188 * If it is positive it must be due to the old information in the
1189 * CG formula, so just remove that and start over with beta=0.
1190 * This corresponds to a steepest descent step.
1195 step--; /* Don't count this step since we are restarting */
1196 continue; /* Go back to the beginning of the big for-loop */
1199 /* Calculate minimum allowed stepsize, before the average (norm)
1200 * relative change in coordinate is smaller than precision
1203 for (i = 0; i < mdatoms->homenr; i++)
1205 for (m = 0; m < DIM; m++)
1207 tmp = fabs(s_min->s.x[i][m]);
1216 /* Add up from all CPUs */
1219 gmx_sumd(1, &minstep, cr);
1222 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1224 if (stepsize < minstep)
1230 /* Write coordinates if necessary */
1231 do_x = do_per_step(step, inputrec->nstxout);
1232 do_f = do_per_step(step, inputrec->nstfout);
1234 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1235 top_global, inputrec, step,
1236 s_min, state_global);
1238 /* Take a step downhill.
1239 * In theory, we should minimize the function along this direction.
1240 * That is quite possible, but it turns out to take 5-10 function evaluations
1241 * for each line. However, we dont really need to find the exact minimum -
1242 * it is much better to start a new CG step in a modified direction as soon
1243 * as we are close to it. This will save a lot of energy evaluations.
1245 * In practice, we just try to take a single step.
1246 * If it worked (i.e. lowered the energy), we increase the stepsize but
1247 * the continue straight to the next CG step without trying to find any minimum.
1248 * If it didn't work (higher energy), there must be a minimum somewhere between
1249 * the old position and the new one.
1251 * Due to the finite numerical accuracy, it turns out that it is a good idea
1252 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1253 * This leads to lower final energies in the tests I've done. / Erik
1255 s_a->epot = s_min->epot;
1257 c = a + stepsize; /* reference position along line is zero */
1259 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1261 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1262 s_min, top, mdatoms, fr, vsite, constr,
1266 /* Take a trial step (new coords in s_c) */
1267 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1268 constr, top, nrnb, wcycle, -1);
1271 /* Calculate energy for the trial step */
1272 evaluate_energy(fplog, cr,
1273 top_global, s_c, top,
1274 inputrec, nrnb, wcycle, gstat,
1275 vsite, constr, fcd, graph, mdatoms, fr,
1276 mu_tot, enerd, vir, pres, -1, FALSE);
1278 /* Calc derivative along line */
1282 for (i = 0; i < mdatoms->homenr; i++)
1284 for (m = 0; m < DIM; m++)
1286 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1289 /* Sum the gradient along the line across CPUs */
1292 gmx_sumd(1, &gpc, cr);
1295 /* This is the max amount of increase in energy we tolerate */
1296 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1298 /* Accept the step if the energy is lower, or if it is not significantly higher
1299 * and the line derivative is still negative.
1301 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1304 /* Great, we found a better energy. Increase step for next iteration
1305 * if we are still going down, decrease it otherwise
1309 stepsize *= 1.618034; /* The golden section */
1313 stepsize *= 0.618034; /* 1/golden section */
1318 /* New energy is the same or higher. We will have to do some work
1319 * to find a smaller value in the interval. Take smaller step next time!
1322 stepsize *= 0.618034;
1328 /* OK, if we didn't find a lower value we will have to locate one now - there must
1329 * be one in the interval [a=0,c].
1330 * The same thing is valid here, though: Don't spend dozens of iterations to find
1331 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1332 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1334 * I also have a safeguard for potentially really pathological functions so we never
1335 * take more than 20 steps before we give up ...
1337 * If we already found a lower value we just skip this step and continue to the update.
1345 /* Select a new trial point.
1346 * If the derivatives at points a & c have different sign we interpolate to zero,
1347 * otherwise just do a bisection.
1349 if (gpa < 0 && gpc > 0)
1351 b = a + gpa*(a-c)/(gpc-gpa);
1358 /* safeguard if interpolation close to machine accuracy causes errors:
1359 * never go outside the interval
1361 if (b <= a || b >= c)
1366 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1368 /* Reload the old state */
1369 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1370 s_min, top, mdatoms, fr, vsite, constr,
1374 /* Take a trial step to this new point - new coords in s_b */
1375 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1376 constr, top, nrnb, wcycle, -1);
1379 /* Calculate energy for the trial step */
1380 evaluate_energy(fplog, cr,
1381 top_global, s_b, top,
1382 inputrec, nrnb, wcycle, gstat,
1383 vsite, constr, fcd, graph, mdatoms, fr,
1384 mu_tot, enerd, vir, pres, -1, FALSE);
1386 /* p does not change within a step, but since the domain decomposition
1387 * might change, we have to use cg_p of s_b here.
1392 for (i = 0; i < mdatoms->homenr; i++)
1394 for (m = 0; m < DIM; m++)
1396 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1399 /* Sum the gradient along the line across CPUs */
1402 gmx_sumd(1, &gpb, cr);
1407 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1408 s_a->epot, s_b->epot, s_c->epot, gpb);
1411 epot_repl = s_b->epot;
1413 /* Keep one of the intervals based on the value of the derivative at the new point */
1416 /* Replace c endpoint with b */
1417 swap_em_state(s_b, s_c);
1423 /* Replace a endpoint with b */
1424 swap_em_state(s_b, s_a);
1430 * Stop search as soon as we find a value smaller than the endpoints.
1431 * Never run more than 20 steps, no matter what.
1435 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1438 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1441 /* OK. We couldn't find a significantly lower energy.
1442 * If beta==0 this was steepest descent, and then we give up.
1443 * If not, set beta=0 and restart with steepest descent before quitting.
1453 /* Reset memory before giving up */
1459 /* Select min energy state of A & C, put the best in B.
1461 if (s_c->epot < s_a->epot)
1465 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1466 s_c->epot, s_a->epot);
1468 swap_em_state(s_b, s_c);
1475 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1476 s_a->epot, s_c->epot);
1478 swap_em_state(s_b, s_a);
1487 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1490 swap_em_state(s_b, s_c);
1494 /* new search direction */
1495 /* beta = 0 means forget all memory and restart with steepest descents. */
1496 if (nstcg && ((step % nstcg) == 0))
1502 /* s_min->fnorm cannot be zero, because then we would have converged
1506 /* Polak-Ribiere update.
1507 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1509 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1511 /* Limit beta to prevent oscillations */
1512 if (fabs(beta) > 5.0)
1518 /* update positions */
1519 swap_em_state(s_min, s_b);
1522 /* Print it if necessary */
1527 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1528 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1529 step, s_min->epot, s_min->fnorm/sqrtNumAtoms,
1530 s_min->fmax, s_min->a_fmax+1);
1533 /* Store the new (lower) energies */
1534 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1535 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1536 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1538 do_log = do_per_step(step, inputrec->nstlog);
1539 do_ene = do_per_step(step, inputrec->nstenergy);
1541 /* Prepare IMD energy record, if bIMD is TRUE. */
1542 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1546 print_ebin_header(fplog, step, step);
1548 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1549 do_log ? fplog : NULL, step, step, eprNORMAL,
1550 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1553 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1554 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1556 IMD_send_positions(inputrec->imd);
1559 /* Stop when the maximum force lies below tolerance.
1560 * If we have reached machine precision, converged is already set to true.
1562 converged = converged || (s_min->fmax < inputrec->em_tol);
1564 } /* End of the loop */
1566 /* IMD cleanup, if bIMD is TRUE. */
1567 IMD_finalize(inputrec->bIMD, inputrec->imd);
1571 step--; /* we never took that last step in this case */
1574 if (s_min->fmax > inputrec->em_tol)
1578 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1579 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1586 /* If we printed energy and/or logfile last step (which was the last step)
1587 * we don't have to do it again, but otherwise print the final values.
1591 /* Write final value to log since we didn't do anything the last step */
1592 print_ebin_header(fplog, step, step);
1594 if (!do_ene || !do_log)
1596 /* Write final energy file entries */
1597 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1598 !do_log ? fplog : NULL, step, step, eprNORMAL,
1599 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1603 /* Print some stuff... */
1606 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1610 * For accurate normal mode calculation it is imperative that we
1611 * store the last conformation into the full precision binary trajectory.
1613 * However, we should only do it if we did NOT already write this step
1614 * above (which we did if do_x or do_f was true).
1616 do_x = !do_per_step(step, inputrec->nstxout);
1617 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1619 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1620 top_global, inputrec, step,
1621 s_min, state_global);
1626 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1627 fnormn = s_min->fnorm/sqrtNumAtoms;
1628 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1629 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1630 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1631 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1633 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1636 finish_em(cr, outf, walltime_accounting, wcycle);
1638 /* To print the actual number of steps we needed somewhere */
1639 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1642 } /* That's all folks */
1645 /*! \brief Do L-BFGS conjugate gradients minimization
1646 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
1647 int nfile, const t_filenm fnm[],
1648 const gmx_output_env_t *oenv, gmx_bool bVerbose,
1650 gmx_vsite_t *vsite, gmx_constr_t constr,
1652 t_inputrec *inputrec,
1653 gmx_mtop_t *top_global, t_fcdata *fcd,
1654 t_state *state_global,
1656 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1659 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
1660 real cpt_period, real max_hours,
1662 unsigned long Flags,
1663 gmx_walltime_accounting_t walltime_accounting)
1665 double do_lbfgs(FILE *fplog, t_commrec *cr,
1666 int nfile, const t_filenm fnm[],
1667 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
1668 int gmx_unused nstglobalcomm,
1669 gmx_vsite_t *vsite, gmx_constr_t constr,
1670 int gmx_unused stepout,
1671 t_inputrec *inputrec,
1672 gmx_mtop_t *top_global, t_fcdata *fcd,
1673 t_state *state_global,
1675 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1676 gmx_edsam_t gmx_unused ed,
1678 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1679 gmx_membed_t gmx_unused *membed,
1680 real gmx_unused cpt_period, real gmx_unused max_hours,
1682 unsigned long gmx_unused Flags,
1683 gmx_walltime_accounting_t walltime_accounting)
1685 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1687 gmx_localtop_t *top;
1688 gmx_enerdata_t *enerd;
1690 gmx_global_stat_t gstat;
1692 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1693 double stepsize, step_taken, gpa, gpb, gpc, tmp, minstep;
1694 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1695 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1696 real a, b, c, maxdelta, delta;
1697 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1698 real dgdx, dgdg, sq, yr, beta;
1703 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1705 int start, end, number_steps;
1707 int i, k, m, n, nfmax, gf, step;
1712 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1717 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).");
1720 n = 3*state_global->natoms;
1721 nmaxcorr = inputrec->nbfgscorr;
1723 /* Allocate memory */
1724 /* Use pointers to real so we dont have to loop over both atoms and
1725 * dimensions all the time...
1726 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1727 * that point to the same memory.
1740 snew(rho, nmaxcorr);
1741 snew(alpha, nmaxcorr);
1744 for (i = 0; i < nmaxcorr; i++)
1750 for (i = 0; i < nmaxcorr; i++)
1759 init_em(fplog, LBFGS, cr, inputrec,
1760 state_global, top_global, &ems, &top, &f,
1761 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1762 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
1763 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1764 * so we free some memory again.
1769 xx = (real *)state_global->x;
1773 end = mdatoms->homenr;
1775 /* Print to log file */
1776 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1778 do_log = do_ene = do_x = do_f = TRUE;
1780 /* Max number of steps */
1781 number_steps = inputrec->nsteps;
1783 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1785 for (i = start; i < end; i++)
1787 if (mdatoms->cFREEZE)
1789 gf = mdatoms->cFREEZE[i];
1791 for (m = 0; m < DIM; m++)
1793 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1798 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1802 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1807 construct_vsites(vsite, state_global->x, 1, NULL,
1808 top->idef.iparams, top->idef.il,
1809 fr->ePBC, fr->bMolPBC, cr, state_global->box);
1812 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1813 /* do_force always puts the charge groups in the box and shifts again
1814 * We do not unshift, so molecules are always whole
1817 ems.s.x = state_global->x;
1819 evaluate_energy(fplog, cr,
1820 top_global, &ems, top,
1821 inputrec, nrnb, wcycle, gstat,
1822 vsite, constr, fcd, graph, mdatoms, fr,
1823 mu_tot, enerd, vir, pres, -1, TRUE);
1828 /* Copy stuff to the energy bin for easy printing etc. */
1829 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1830 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
1831 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1833 print_ebin_header(fplog, step, step);
1834 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1835 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1839 /* This is the starting energy */
1840 Epot = enerd->term[F_EPOT];
1846 /* Set the initial step.
1847 * since it will be multiplied by the non-normalized search direction
1848 * vector (force vector the first time), we scale it by the
1849 * norm of the force.
1854 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
1855 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1856 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1857 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrtNumAtoms);
1858 fprintf(stderr, "\n");
1859 /* and copy to the log file too... */
1860 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1861 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1862 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrtNumAtoms);
1863 fprintf(fplog, "\n");
1866 // Point is an index to the memory of search directions, where 0 is the first one.
1869 // Set initial search direction to the force (-gradient), or 0 for frozen particles.
1870 for (i = 0; i < n; i++)
1874 dx[point][i] = ff[i]; /* Initial search direction */
1882 // Stepsize will be modified during the search, and actually it is not critical
1883 // (the main efficiency in the algorithm comes from changing directions), but
1884 // we still need an initial value, so estimate it as the inverse of the norm
1885 // so we take small steps where the potential fluctuates a lot.
1886 stepsize = 1.0/fnorm;
1888 /* Start the loop over BFGS steps.
1889 * Each successful step is counted, and we continue until
1890 * we either converge or reach the max number of steps.
1895 /* Set the gradient from the force */
1897 for (step = 0; (number_steps < 0 || step <= number_steps) && !converged; step++)
1900 /* Write coordinates if necessary */
1901 do_x = do_per_step(step, inputrec->nstxout);
1902 do_f = do_per_step(step, inputrec->nstfout);
1907 mdof_flags |= MDOF_X;
1912 mdof_flags |= MDOF_F;
1917 mdof_flags |= MDOF_IMD;
1920 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1921 top_global, step, (real)step, state_global, state_global, f);
1923 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1925 /* make s a pointer to current search direction - point=0 first time we get here */
1928 // calculate line gradient in position A
1929 for (gpa = 0, i = 0; i < n; i++)
1934 /* Calculate minimum allowed stepsize along the line, before the average (norm)
1935 * relative change in coordinate is smaller than precision
1937 for (minstep = 0, i = 0; i < n; i++)
1947 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1949 if (stepsize < minstep)
1955 // Before taking any steps along the line, store the old position
1956 for (i = 0; i < n; i++)
1963 for (i = 0; i < n; i++)
1968 /* Take a step downhill.
1969 * In theory, we should find the actual minimum of the function in this
1970 * direction, somewhere along the line.
1971 * That is quite possible, but it turns out to take 5-10 function evaluations
1972 * for each line. However, we dont really need to find the exact minimum -
1973 * it is much better to start a new BFGS step in a modified direction as soon
1974 * as we are close to it. This will save a lot of energy evaluations.
1976 * In practice, we just try to take a single step.
1977 * If it worked (i.e. lowered the energy), we increase the stepsize but
1978 * continue straight to the next BFGS step without trying to find any minimum,
1979 * i.e. we change the search direction too. If the line was smooth, it is
1980 * likely we are in a smooth region, and then it makes sense to take longer
1981 * steps in the modified search direction too.
1983 * If it didn't work (higher energy), there must be a minimum somewhere between
1984 * the old position and the new one. Then we need to start by finding a lower
1985 * value before we change search direction. Since the energy was apparently
1986 * quite rough, we need to decrease the step size.
1988 * Due to the finite numerical accuracy, it turns out that it is a good idea
1989 * to accept a SMALL increase in energy, if the derivative is still downhill.
1990 * This leads to lower final energies in the tests I've done. / Erik
1993 // State "A" is the first position along the line.
1994 // reference position along line is initially zero
1998 // Check stepsize first. We do not allow displacements
1999 // larger than emstep.
2003 // Pick a new position C by adding stepsize to A.
2006 // Calculate what the largest change in any individual coordinate
2007 // would be (translation along line * gradient along line)
2009 for (i = 0; i < n; i++)
2012 if (delta > maxdelta)
2017 // If any displacement is larger than the stepsize limit, reduce the step
2018 if (maxdelta > inputrec->em_stepsize)
2023 while (maxdelta > inputrec->em_stepsize);
2025 // Take a trial step and move the coordinate array xc[] to position C
2026 for (i = 0; i < n; i++)
2028 xc[i] = lastx[i] + c*s[i];
2032 // Calculate energy for the trial step in position C
2033 ems.s.x = (rvec *)xc;
2035 evaluate_energy(fplog, cr,
2036 top_global, &ems, top,
2037 inputrec, nrnb, wcycle, gstat,
2038 vsite, constr, fcd, graph, mdatoms, fr,
2039 mu_tot, enerd, vir, pres, step, FALSE);
2042 // Calc line gradient in position C
2043 for (gpc = 0, i = 0; i < n; i++)
2045 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
2047 /* Sum the gradient along the line across CPUs */
2050 gmx_sumd(1, &gpc, cr);
2053 // This is the max amount of increase in energy we tolerate.
2054 // By allowing VERY small changes (close to numerical precision) we
2055 // frequently find even better (lower) final energies.
2056 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
2058 // Accept the step if the energy is lower in the new position C (compared to A),
2059 // or if it is not significantly higher and the line derivative is still negative.
2060 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
2062 // Great, we found a better energy. We no longer try to alter the
2063 // stepsize, but simply accept this new better position. The we select a new
2064 // search direction instead, which will be much more efficient than continuing
2065 // to take smaller steps along a line. Set fnorm based on the new C position,
2066 // which will be used to update the stepsize to 1/fnorm further down.
2072 // If we got here, the energy is NOT lower in point C, i.e. it will be the same
2073 // or higher than in point A. In this case it is pointless to move to point C,
2074 // so we will have to do more iterations along the same line to find a smaller
2075 // value in the interval [A=0.0,C].
2076 // Here, A is still 0.0, but that will change when we do a search in the interval
2077 // [0.0,C] below. That search we will do by interpolation or bisection rather
2078 // than with the stepsize, so no need to modify it. For the next search direction
2079 // it will be reset to 1/fnorm anyway.
2085 // OK, if we didn't find a lower value we will have to locate one now - there must
2086 // be one in the interval [a,c].
2087 // The same thing is valid here, though: Don't spend dozens of iterations to find
2088 // the line minimum. We try to interpolate based on the derivative at the endpoints,
2089 // and only continue until we find a lower value. In most cases this means 1-2 iterations.
2090 // I also have a safeguard for potentially really pathological functions so we never
2091 // take more than 20 steps before we give up.
2092 // If we already found a lower value we just skip this step and continue to the update.
2096 // Select a new trial point B in the interval [A,C].
2097 // If the derivatives at points a & c have different sign we interpolate to zero,
2098 // otherwise just do a bisection since there might be multiple minima/maxima
2099 // inside the interval.
2100 if (gpa < 0 && gpc > 0)
2102 b = a + gpa*(a-c)/(gpc-gpa);
2109 /* safeguard if interpolation close to machine accuracy causes errors:
2110 * never go outside the interval
2112 if (b <= a || b >= c)
2117 // Take a trial step to point B
2118 for (i = 0; i < n; i++)
2120 xb[i] = lastx[i] + b*s[i];
2124 // Calculate energy for the trial step in point B
2125 ems.s.x = (rvec *)xb;
2127 evaluate_energy(fplog, cr,
2128 top_global, &ems, top,
2129 inputrec, nrnb, wcycle, gstat,
2130 vsite, constr, fcd, graph, mdatoms, fr,
2131 mu_tot, enerd, vir, pres, step, FALSE);
2135 // Calculate gradient in point B
2136 for (gpb = 0, i = 0; i < n; i++)
2138 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2141 /* Sum the gradient along the line across CPUs */
2144 gmx_sumd(1, &gpb, cr);
2147 // Keep one of the intervals [A,B] or [B,C] based on the value of the derivative
2148 // at the new point B, and rename the endpoints of this new interval A and C.
2151 /* Replace c endpoint with b */
2155 /* swap coord pointers b/c */
2165 /* Replace a endpoint with b */
2169 /* swap coord pointers a/b */
2179 * Stop search as soon as we find a value smaller than the endpoints,
2180 * or if the tolerance is below machine precision.
2181 * Never run more than 20 steps, no matter what.
2185 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2187 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2189 /* OK. We couldn't find a significantly lower energy.
2190 * If ncorr==0 this was steepest descent, and then we give up.
2191 * If not, reset memory to restart as steepest descent before quitting.
2203 /* Search in gradient direction */
2204 for (i = 0; i < n; i++)
2206 dx[point][i] = ff[i];
2208 /* Reset stepsize */
2209 stepsize = 1.0/fnorm;
2214 /* Select min energy state of A & C, put the best in xx/ff/Epot
2220 for (i = 0; i < n; i++)
2231 for (i = 0; i < n; i++)
2245 for (i = 0; i < n; i++)
2253 /* Update the memory information, and calculate a new
2254 * approximation of the inverse hessian
2257 /* Have new data in Epot, xx, ff */
2258 if (ncorr < nmaxcorr)
2263 for (i = 0; i < n; i++)
2265 dg[point][i] = lastf[i]-ff[i];
2266 dx[point][i] *= step_taken;
2271 for (i = 0; i < n; i++)
2273 dgdg += dg[point][i]*dg[point][i];
2274 dgdx += dg[point][i]*dx[point][i];
2279 rho[point] = 1.0/dgdx;
2282 if (point >= nmaxcorr)
2288 for (i = 0; i < n; i++)
2295 /* Recursive update. First go back over the memory points */
2296 for (k = 0; k < ncorr; k++)
2305 for (i = 0; i < n; i++)
2307 sq += dx[cp][i]*p[i];
2310 alpha[cp] = rho[cp]*sq;
2312 for (i = 0; i < n; i++)
2314 p[i] -= alpha[cp]*dg[cp][i];
2318 for (i = 0; i < n; i++)
2323 /* And then go forward again */
2324 for (k = 0; k < ncorr; k++)
2327 for (i = 0; i < n; i++)
2329 yr += p[i]*dg[cp][i];
2333 beta = alpha[cp]-beta;
2335 for (i = 0; i < n; i++)
2337 p[i] += beta*dx[cp][i];
2347 for (i = 0; i < n; i++)
2351 dx[point][i] = p[i];
2359 /* Test whether the convergence criterion is met */
2360 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2362 /* Print it if necessary */
2367 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2368 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2369 step, Epot, fnorm/sqrtNumAtoms, fmax, nfmax+1);
2372 /* Store the new (lower) energies */
2373 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2374 mdatoms->tmass, enerd, state_global, inputrec->fepvals, inputrec->expandedvals, state_global->box,
2375 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2376 do_log = do_per_step(step, inputrec->nstlog);
2377 do_ene = do_per_step(step, inputrec->nstenergy);
2380 print_ebin_header(fplog, step, step);
2382 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2383 do_log ? fplog : NULL, step, step, eprNORMAL,
2384 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2387 /* Send x and E to IMD client, if bIMD is TRUE. */
2388 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2390 IMD_send_positions(inputrec->imd);
2393 // Reset stepsize in we are doing more iterations
2394 stepsize = 1.0/fnorm;
2396 /* Stop when the maximum force lies below tolerance.
2397 * If we have reached machine precision, converged is already set to true.
2399 converged = converged || (fmax < inputrec->em_tol);
2401 } /* End of the loop */
2403 /* IMD cleanup, if bIMD is TRUE. */
2404 IMD_finalize(inputrec->bIMD, inputrec->imd);
2408 step--; /* we never took that last step in this case */
2411 if (fmax > inputrec->em_tol)
2415 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2416 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2421 /* If we printed energy and/or logfile last step (which was the last step)
2422 * we don't have to do it again, but otherwise print the final values.
2424 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2426 print_ebin_header(fplog, step, step);
2428 if (!do_ene || !do_log) /* Write final energy file entries */
2430 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2431 !do_log ? fplog : NULL, step, step, eprNORMAL,
2432 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2435 /* Print some stuff... */
2438 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2442 * For accurate normal mode calculation it is imperative that we
2443 * store the last conformation into the full precision binary trajectory.
2445 * However, we should only do it if we did NOT already write this step
2446 * above (which we did if do_x or do_f was true).
2448 do_x = !do_per_step(step, inputrec->nstxout);
2449 do_f = !do_per_step(step, inputrec->nstfout);
2450 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2451 top_global, inputrec, step,
2452 &ems, state_global);
2456 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2457 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2458 number_steps, Epot, fmax, nfmax, fnorm/sqrtNumAtoms);
2459 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2460 number_steps, Epot, fmax, nfmax, fnorm/sqrtNumAtoms);
2462 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2465 finish_em(cr, outf, walltime_accounting, wcycle);
2467 /* To print the actual number of steps we needed somewhere */
2468 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2471 } /* That's all folks */
2473 /*! \brief Do steepest descents minimization
2474 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
2475 int nfile, const t_filenm fnm[],
2476 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2478 gmx_vsite_t *vsite, gmx_constr_t constr,
2480 t_inputrec *inputrec,
2481 gmx_mtop_t *top_global, t_fcdata *fcd,
2482 t_state *state_global,
2484 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2487 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2488 real cpt_period, real max_hours,
2490 unsigned long Flags,
2491 gmx_walltime_accounting_t walltime_accounting)
2493 double do_steep(FILE *fplog, t_commrec *cr,
2494 int nfile, const t_filenm fnm[],
2495 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2496 int gmx_unused nstglobalcomm,
2497 gmx_vsite_t *vsite, gmx_constr_t constr,
2498 int gmx_unused stepout,
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,
2504 gmx_edsam_t gmx_unused ed,
2506 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2507 gmx_membed_t gmx_unused *membed,
2508 real gmx_unused cpt_period, real gmx_unused max_hours,
2510 unsigned long gmx_unused Flags,
2511 gmx_walltime_accounting_t walltime_accounting)
2513 const char *SD = "Steepest Descents";
2514 em_state_t *s_min, *s_try;
2515 gmx_localtop_t *top;
2516 gmx_enerdata_t *enerd;
2518 gmx_global_stat_t gstat;
2524 gmx_bool bDone, bAbort, do_x, do_f;
2529 int steps_accepted = 0;
2531 s_min = init_em_state();
2532 s_try = init_em_state();
2534 /* Init em and store the local state in s_try */
2535 init_em(fplog, SD, cr, inputrec,
2536 state_global, top_global, s_try, &top, &f,
2537 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2538 nfile, fnm, &outf, &mdebin, imdport, Flags, wcycle);
2540 /* Print to log file */
2541 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2543 /* Set variables for stepsize (in nm). This is the largest
2544 * step that we are going to make in any direction.
2546 ustep = inputrec->em_stepsize;
2549 /* Max number of steps */
2550 nsteps = inputrec->nsteps;
2554 /* Print to the screen */
2555 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2559 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2562 /**** HERE STARTS THE LOOP ****
2563 * count is the counter for the number of steps
2564 * bDone will be TRUE when the minimization has converged
2565 * bAbort will be TRUE when nsteps steps have been performed or when
2566 * the stepsize becomes smaller than is reasonable for machine precision
2571 while (!bDone && !bAbort)
2573 bAbort = (nsteps >= 0) && (count == nsteps);
2575 /* set new coordinates, except for first step */
2576 bool validStep = true;
2580 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2581 s_min, stepsize, s_min->f, s_try,
2582 constr, top, nrnb, wcycle, count);
2587 evaluate_energy(fplog, cr,
2588 top_global, s_try, top,
2589 inputrec, nrnb, wcycle, gstat,
2590 vsite, constr, fcd, graph, mdatoms, fr,
2591 mu_tot, enerd, vir, pres, count, count == 0);
2595 // Signal constraint error during stepping with energy=inf
2596 s_try->epot = std::numeric_limits<real>::infinity();
2601 print_ebin_header(fplog, count, count);
2606 s_min->epot = s_try->epot;
2609 /* Print it if necessary */
2614 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2615 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2616 ( (count == 0) || (s_try->epot < s_min->epot) ) ? '\n' : '\r');
2620 if ( (count == 0) || (s_try->epot < s_min->epot) )
2622 /* Store the new (lower) energies */
2623 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2624 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2625 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2627 /* Prepare IMD energy record, if bIMD is TRUE. */
2628 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2630 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2631 do_per_step(steps_accepted, inputrec->nstdisreout),
2632 do_per_step(steps_accepted, inputrec->nstorireout),
2633 fplog, count, count, eprNORMAL,
2634 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2639 /* Now if the new energy is smaller than the previous...
2640 * or if this is the first step!
2641 * or if we did random steps!
2644 if ( (count == 0) || (s_try->epot < s_min->epot) )
2648 /* Test whether the convergence criterion is met... */
2649 bDone = (s_try->fmax < inputrec->em_tol);
2651 /* Copy the arrays for force, positions and energy */
2652 /* The 'Min' array always holds the coords and forces of the minimal
2654 swap_em_state(s_min, s_try);
2660 /* Write to trn, if necessary */
2661 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2662 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2663 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2664 top_global, inputrec, count,
2665 s_min, state_global);
2669 /* If energy is not smaller make the step smaller... */
2672 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2674 /* Reload the old state */
2675 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2676 s_min, top, mdatoms, fr, vsite, constr,
2681 /* Determine new step */
2682 stepsize = ustep/s_min->fmax;
2684 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2686 if (count == nsteps || ustep < 1e-12)
2688 if (count == nsteps || ustep < 1e-6)
2693 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2694 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2699 /* Send IMD energies and positions, if bIMD is TRUE. */
2700 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2702 IMD_send_positions(inputrec->imd);
2706 } /* End of the loop */
2708 /* IMD cleanup, if bIMD is TRUE. */
2709 IMD_finalize(inputrec->bIMD, inputrec->imd);
2711 /* Print some data... */
2714 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2716 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2717 top_global, inputrec, count,
2718 s_min, state_global);
2722 double sqrtNumAtoms = sqrt(static_cast<double>(state_global->natoms));
2723 fnormn = s_min->fnorm/sqrtNumAtoms;
2725 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2726 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2727 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2728 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2731 finish_em(cr, outf, walltime_accounting, wcycle);
2733 /* To print the actual number of steps we needed somewhere */
2734 inputrec->nsteps = count;
2736 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2739 } /* That's all folks */
2741 /*! \brief Do normal modes analysis
2742 \copydoc integrator_t (FILE *fplog, t_commrec *cr,
2743 int nfile, const t_filenm fnm[],
2744 const gmx_output_env_t *oenv, gmx_bool bVerbose,
2746 gmx_vsite_t *vsite, gmx_constr_t constr,
2748 t_inputrec *inputrec,
2749 gmx_mtop_t *top_global, t_fcdata *fcd,
2750 t_state *state_global,
2752 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2755 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
2756 real cpt_period, real max_hours,
2758 unsigned long Flags,
2759 gmx_walltime_accounting_t walltime_accounting)
2761 double do_nm(FILE *fplog, t_commrec *cr,
2762 int nfile, const t_filenm fnm[],
2763 const gmx_output_env_t gmx_unused *oenv, gmx_bool bVerbose,
2764 int gmx_unused nstglobalcomm,
2765 gmx_vsite_t *vsite, gmx_constr_t constr,
2766 int gmx_unused stepout,
2767 t_inputrec *inputrec,
2768 gmx_mtop_t *top_global, t_fcdata *fcd,
2769 t_state *state_global,
2771 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2772 gmx_edsam_t gmx_unused ed,
2774 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2775 gmx_membed_t gmx_unused *membed,
2776 real gmx_unused cpt_period, real gmx_unused max_hours,
2778 unsigned long gmx_unused Flags,
2779 gmx_walltime_accounting_t walltime_accounting)
2781 const char *NM = "Normal Mode Analysis";
2784 gmx_localtop_t *top;
2785 gmx_enerdata_t *enerd;
2787 gmx_global_stat_t gstat;
2792 gmx_bool bSparse; /* use sparse matrix storage format */
2794 gmx_sparsematrix_t * sparse_matrix = NULL;
2795 real * full_matrix = NULL;
2796 em_state_t * state_work;
2798 /* added with respect to mdrun */
2800 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2802 bool bIsMaster = MASTER(cr);
2806 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2809 state_work = init_em_state();
2811 /* Init em and store the local state in state_minimum */
2812 init_em(fplog, NM, cr, inputrec,
2813 state_global, top_global, state_work, &top,
2815 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2816 nfile, fnm, &outf, NULL, imdport, Flags, wcycle);
2818 gmx_shellfc_t *shellfc = init_shell_flexcon(stdout,
2820 n_flexible_constraints(constr),
2821 inputrec->nstcalcenergy,
2826 make_local_shells(cr, mdatoms, shellfc);
2828 std::vector<size_t> atom_index = get_atom_index(top_global);
2829 snew(fneg, atom_index.size());
2830 snew(dfdx, atom_index.size());
2836 "NOTE: This version of GROMACS has been compiled in single precision,\n"
2837 " which MIGHT not be accurate enough for normal mode analysis.\n"
2838 " GROMACS now uses sparse matrix storage, so the memory requirements\n"
2839 " are fairly modest even if you recompile in double precision.\n\n");
2843 /* Check if we can/should use sparse storage format.
2845 * Sparse format is only useful when the Hessian itself is sparse, which it
2846 * will be when we use a cutoff.
2847 * For small systems (n<1000) it is easier to always use full matrix format, though.
2849 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2851 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2854 else if (atom_index.size() < 1000)
2856 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", atom_index.size());
2861 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2865 /* Number of dimensions, based on real atoms, that is not vsites or shell */
2866 sz = DIM*atom_index.size();
2868 fprintf(stderr, "Allocating Hessian memory...\n\n");
2872 sparse_matrix = gmx_sparsematrix_init(sz);
2873 sparse_matrix->compressed_symmetric = TRUE;
2877 snew(full_matrix, sz*sz);
2884 /* Write start time and temperature */
2885 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2887 /* fudge nr of steps to nr of atoms */
2888 inputrec->nsteps = atom_index.size()*2;
2892 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2893 *(top_global->name), (int)inputrec->nsteps);
2896 nnodes = cr->nnodes;
2898 /* Make evaluate_energy do a single node force calculation */
2900 evaluate_energy(fplog, cr,
2901 top_global, state_work, top,
2902 inputrec, nrnb, wcycle, gstat,
2903 vsite, constr, fcd, graph, mdatoms, fr,
2904 mu_tot, enerd, vir, pres, -1, TRUE);
2905 cr->nnodes = nnodes;
2907 /* if forces are not small, warn user */
2908 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2910 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2911 if (state_work->fmax > 1.0e-3)
2913 md_print_info(cr, fplog,
2914 "The force is probably not small enough to "
2915 "ensure that you are at a minimum.\n"
2916 "Be aware that negative eigenvalues may occur\n"
2917 "when the resulting matrix is diagonalized.\n\n");
2920 /***********************************************************
2922 * Loop over all pairs in matrix
2924 * do_force called twice. Once with positive and
2925 * once with negative displacement
2927 ************************************************************/
2929 /* Steps are divided one by one over the nodes */
2931 for (unsigned int aid = cr->nodeid; aid < atom_index.size(); aid += nnodes)
2933 size_t atom = atom_index[aid];
2934 for (size_t d = 0; d < DIM; d++)
2936 gmx_bool bBornRadii = FALSE;
2937 gmx_int64_t step = 0;
2938 int force_flags = GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES;
2941 x_min = state_work->s.x[atom][d];
2943 for (unsigned int dx = 0; (dx < 2); dx++)
2947 state_work->s.x[atom][d] = x_min - der_range;
2951 state_work->s.x[atom][d] = x_min + der_range;
2954 /* Make evaluate_energy do a single node force calculation */
2958 /* Now is the time to relax the shells */
2959 (void) relax_shell_flexcon(fplog, cr, bVerbose, step,
2960 inputrec, bNS, force_flags,
2963 &state_work->s, state_work->f, vir, mdatoms,
2964 nrnb, wcycle, graph, &top_global->groups,
2965 shellfc, fr, bBornRadii, t, mu_tot,
2972 evaluate_energy(fplog, cr,
2973 top_global, state_work, top,
2974 inputrec, nrnb, wcycle, gstat,
2975 vsite, constr, fcd, graph, mdatoms, fr,
2976 mu_tot, enerd, vir, pres, atom*2+dx, FALSE);
2979 cr->nnodes = nnodes;
2983 for (size_t i = 0; i < atom_index.size(); i++)
2985 copy_rvec(state_work->f[atom_index[i]], fneg[i]);
2990 /* x is restored to original */
2991 state_work->s.x[atom][d] = x_min;
2993 for (size_t j = 0; j < atom_index.size(); j++)
2995 for (size_t k = 0; (k < DIM); k++)
2998 -(state_work->f[atom_index[j]][k] - fneg[j][k])/(2*der_range);
3005 #define mpi_type GMX_MPI_REAL
3006 MPI_Send(dfdx[0], atom_index.size()*DIM, mpi_type, MASTER(cr),
3007 cr->nodeid, cr->mpi_comm_mygroup);
3012 for (node = 0; (node < nnodes && atom+node < atom_index.size()); node++)
3018 MPI_Recv(dfdx[0], atom_index.size()*DIM, mpi_type, node, node,
3019 cr->mpi_comm_mygroup, &stat);
3024 row = (atom + node)*DIM + d;
3026 for (size_t j = 0; j < atom_index.size(); j++)
3028 for (size_t k = 0; k < DIM; k++)
3034 if (col >= row && dfdx[j][k] != 0.0)
3036 gmx_sparsematrix_increment_value(sparse_matrix,
3037 row, col, dfdx[j][k]);
3042 full_matrix[row*sz+col] = dfdx[j][k];
3049 if (bVerbose && fplog)
3054 /* write progress */
3055 if (bIsMaster && bVerbose)
3057 fprintf(stderr, "\rFinished step %d out of %d",
3058 static_cast<int>(std::min(atom+nnodes, atom_index.size())),
3059 static_cast<int>(atom_index.size()));
3066 fprintf(stderr, "\n\nWriting Hessian...\n");
3067 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
3070 finish_em(cr, outf, walltime_accounting, wcycle);
3072 walltime_accounting_set_nsteps_done(walltime_accounting, atom_index.size()*2);