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, 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.
57 #include "md_support.h"
62 #include "gromacs/topology/mtop_util.h"
65 #include "gmx_omp_nthreads.h"
66 #include "md_logging.h"
68 #include "gromacs/fileio/confio.h"
69 #include "gromacs/fileio/mtxio.h"
70 #include "gromacs/fileio/trajectory_writing.h"
71 #include "gromacs/imd/imd.h"
72 #include "gromacs/legacyheaders/types/commrec.h"
73 #include "gromacs/linearalgebra/sparsematrix.h"
74 #include "gromacs/math/vec.h"
75 #include "gromacs/pbcutil/mshift.h"
76 #include "gromacs/pbcutil/pbc.h"
77 #include "gromacs/timing/wallcycle.h"
78 #include "gromacs/timing/walltime_accounting.h"
79 #include "gromacs/utility/cstringutil.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/smalloc.h"
92 static em_state_t *init_em_state()
98 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
99 snew(ems->s.lambda, efptNR);
104 static void print_em_start(FILE *fplog,
106 gmx_walltime_accounting_t walltime_accounting,
107 gmx_wallcycle_t wcycle,
110 walltime_accounting_start(walltime_accounting);
111 wallcycle_start(wcycle, ewcRUN);
112 print_start(fplog, cr, walltime_accounting, name);
114 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
115 gmx_wallcycle_t wcycle)
117 wallcycle_stop(wcycle, ewcRUN);
119 walltime_accounting_end(walltime_accounting);
122 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
125 fprintf(out, "%s:\n", minimizer);
126 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
127 fprintf(out, " Number of steps = %12d\n", nsteps);
130 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
136 "\nEnergy minimization reached the maximum number "
137 "of steps before the forces reached the requested "
138 "precision Fmax < %g.\n", ftol);
143 "\nEnergy minimization has stopped, but the forces have "
144 "not converged to the requested precision Fmax < %g (which "
145 "may not be possible for your system). It stopped "
146 "because the algorithm tried to make a new step whose size "
147 "was too small, or there was no change in the energy since "
148 "last step. Either way, we regard the minimization as "
149 "converged to within the available machine precision, "
150 "given your starting configuration and EM parameters.\n%s%s",
152 sizeof(real) < sizeof(double) ?
153 "\nDouble precision normally gives you higher accuracy, but "
154 "this is often not needed for preparing to run molecular "
158 "You might need to increase your constraint accuracy, or turn\n"
159 "off constraints altogether (set constraints = none in mdp file)\n" :
162 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
167 static void print_converged(FILE *fp, const char *alg, real ftol,
168 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
169 real epot, real fmax, int nfmax, real fnorm)
171 char buf[STEPSTRSIZE];
175 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
176 alg, ftol, gmx_step_str(count, buf));
178 else if (count < nsteps)
180 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
181 "but did not reach the requested Fmax < %g.\n",
182 alg, gmx_step_str(count, buf), ftol);
186 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
187 alg, ftol, gmx_step_str(count, buf));
191 fprintf(fp, "Potential Energy = %21.14e\n", epot);
192 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
193 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
195 fprintf(fp, "Potential Energy = %14.7e\n", epot);
196 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
197 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
201 static void get_f_norm_max(t_commrec *cr,
202 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
203 real *fnorm, real *fmax, int *a_fmax)
206 real fmax2, fmax2_0, fam;
207 int la_max, a_max, start, end, i, m, gf;
209 /* This routine finds the largest force and returns it.
210 * On parallel machines the global max is taken.
217 end = mdatoms->homenr;
218 if (mdatoms->cFREEZE)
220 for (i = start; i < end; i++)
222 gf = mdatoms->cFREEZE[i];
224 for (m = 0; m < DIM; m++)
226 if (!opts->nFreeze[gf][m])
241 for (i = start; i < end; i++)
253 if (la_max >= 0 && DOMAINDECOMP(cr))
255 a_max = cr->dd->gatindex[la_max];
263 snew(sum, 2*cr->nnodes+1);
264 sum[2*cr->nodeid] = fmax2;
265 sum[2*cr->nodeid+1] = a_max;
266 sum[2*cr->nnodes] = fnorm2;
267 gmx_sumd(2*cr->nnodes+1, sum, cr);
268 fnorm2 = sum[2*cr->nnodes];
269 /* Determine the global maximum */
270 for (i = 0; i < cr->nnodes; i++)
272 if (sum[2*i] > fmax2)
275 a_max = (int)(sum[2*i+1] + 0.5);
283 *fnorm = sqrt(fnorm2);
295 static void get_state_f_norm_max(t_commrec *cr,
296 t_grpopts *opts, t_mdatoms *mdatoms,
299 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
302 void init_em(FILE *fplog, const char *title,
303 t_commrec *cr, t_inputrec *ir,
304 t_state *state_global, gmx_mtop_t *top_global,
305 em_state_t *ems, gmx_localtop_t **top,
306 rvec **f, rvec **f_global,
307 t_nrnb *nrnb, rvec mu_tot,
308 t_forcerec *fr, gmx_enerdata_t **enerd,
309 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
310 gmx_vsite_t *vsite, gmx_constr_t constr,
311 int nfile, const t_filenm fnm[],
312 gmx_mdoutf_t *outf, t_mdebin **mdebin,
313 int imdport, unsigned long gmx_unused Flags)
320 fprintf(fplog, "Initiating %s\n", title);
323 state_global->ngtc = 0;
325 /* Initialize lambda variables */
326 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
330 /* Interactive molecular dynamics */
331 init_IMD(ir, cr, top_global, fplog, 1, state_global->x,
332 nfile, fnm, NULL, imdport, Flags);
334 if (DOMAINDECOMP(cr))
336 *top = dd_init_local_top(top_global);
338 dd_init_local_state(cr->dd, state_global, &ems->s);
342 /* Distribute the charge groups over the nodes from the master node */
343 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
344 state_global, top_global, ir,
345 &ems->s, &ems->f, mdatoms, *top,
346 fr, vsite, NULL, constr,
348 dd_store_state(cr->dd, &ems->s);
352 snew(*f_global, top_global->natoms);
362 snew(*f, top_global->natoms);
364 /* Just copy the state */
365 ems->s = *state_global;
366 snew(ems->s.x, ems->s.nalloc);
367 snew(ems->f, ems->s.nalloc);
368 for (i = 0; i < state_global->natoms; i++)
370 copy_rvec(state_global->x[i], ems->s.x[i]);
372 copy_mat(state_global->box, ems->s.box);
374 *top = gmx_mtop_generate_local_top(top_global, ir);
377 forcerec_set_excl_load(fr, *top);
379 setup_bonded_threading(fr, &(*top)->idef);
381 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
383 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
390 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
391 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
395 set_vsite_top(vsite, *top, mdatoms, cr);
401 if (ir->eConstrAlg == econtSHAKE &&
402 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
404 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
405 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
408 if (!DOMAINDECOMP(cr))
410 set_constraints(constr, *top, ir, mdatoms, cr);
413 if (!ir->bContinuation)
415 /* Constrain the starting coordinates */
417 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
418 ir, NULL, cr, -1, 0, 1.0, mdatoms,
419 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
420 ems->s.lambda[efptFEP], &dvdl_constr,
421 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
427 *gstat = global_stat_init(ir);
430 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, NULL);
433 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
438 /* Init bin for energy stuff */
439 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
443 calc_shifts(ems->s.box, fr->shift_vec);
446 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
447 gmx_walltime_accounting_t walltime_accounting,
448 gmx_wallcycle_t wcycle)
450 if (!(cr->duty & DUTY_PME))
452 /* Tell the PME only node to finish */
453 gmx_pme_send_finish(cr);
458 em_time_end(walltime_accounting, wcycle);
461 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
470 static void copy_em_coords(em_state_t *ems, t_state *state)
474 for (i = 0; (i < state->natoms); i++)
476 copy_rvec(ems->s.x[i], state->x[i]);
480 static void write_em_traj(FILE *fplog, t_commrec *cr,
482 gmx_bool bX, gmx_bool bF, const char *confout,
483 gmx_mtop_t *top_global,
484 t_inputrec *ir, gmx_int64_t step,
486 t_state *state_global, rvec *f_global)
489 gmx_bool bIMDout = FALSE;
492 /* Shall we do IMD output? */
495 bIMDout = do_per_step(step, IMD_get_step(ir->imd->setup));
498 if ((bX || bF || bIMDout || confout != NULL) && !DOMAINDECOMP(cr))
500 copy_em_coords(state, state_global);
507 mdof_flags |= MDOF_X;
511 mdof_flags |= MDOF_F;
514 /* If we want IMD output, set appropriate MDOF flag */
517 mdof_flags |= MDOF_IMD;
520 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
521 top_global, step, (double)step,
522 &state->s, state_global, state->f, f_global);
524 if (confout != NULL && MASTER(cr))
526 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
528 /* Make molecules whole only for confout writing */
529 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
533 write_sto_conf_mtop(confout,
534 *top_global->name, top_global,
535 state_global->x, NULL, ir->ePBC, state_global->box);
539 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
541 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
542 gmx_constr_t constr, gmx_localtop_t *top,
543 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
556 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
558 gmx_incons("state mismatch in do_em_step");
561 s2->flags = s1->flags;
563 if (s2->nalloc != s1->nalloc)
565 s2->nalloc = s1->nalloc;
566 srenew(s2->x, s1->nalloc);
567 srenew(ems2->f, s1->nalloc);
568 if (s2->flags & (1<<estCGP))
570 srenew(s2->cg_p, s1->nalloc);
574 s2->natoms = s1->natoms;
575 copy_mat(s1->box, s2->box);
576 /* Copy free energy state */
577 for (i = 0; i < efptNR; i++)
579 s2->lambda[i] = s1->lambda[i];
581 copy_mat(s1->box, s2->box);
589 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
594 #pragma omp for schedule(static) nowait
595 for (i = start; i < end; i++)
601 for (m = 0; m < DIM; m++)
603 if (ir->opts.nFreeze[gf][m])
609 x2[i][m] = x1[i][m] + a*f[i][m];
614 if (s2->flags & (1<<estCGP))
616 /* Copy the CG p vector */
619 #pragma omp for schedule(static) nowait
620 for (i = start; i < end; i++)
622 copy_rvec(x1[i], x2[i]);
626 if (DOMAINDECOMP(cr))
628 s2->ddp_count = s1->ddp_count;
629 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
632 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
633 srenew(s2->cg_gl, s2->cg_gl_nalloc);
636 s2->ncg_gl = s1->ncg_gl;
637 #pragma omp for schedule(static) nowait
638 for (i = 0; i < s2->ncg_gl; i++)
640 s2->cg_gl[i] = s1->cg_gl[i];
642 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
648 wallcycle_start(wcycle, ewcCONSTR);
650 constrain(NULL, TRUE, TRUE, constr, &top->idef,
651 ir, NULL, cr, count, 0, 1.0, md,
652 s1->x, s2->x, NULL, bMolPBC, s2->box,
653 s2->lambda[efptBONDED], &dvdl_constr,
654 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
655 wallcycle_stop(wcycle, ewcCONSTR);
659 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
660 gmx_mtop_t *top_global, t_inputrec *ir,
661 em_state_t *ems, gmx_localtop_t *top,
662 t_mdatoms *mdatoms, t_forcerec *fr,
663 gmx_vsite_t *vsite, gmx_constr_t constr,
664 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
666 /* Repartition the domain decomposition */
667 wallcycle_start(wcycle, ewcDOMDEC);
668 dd_partition_system(fplog, step, cr, FALSE, 1,
669 NULL, top_global, ir,
671 mdatoms, top, fr, vsite, NULL, constr,
672 nrnb, wcycle, FALSE);
673 dd_store_state(cr->dd, &ems->s);
674 wallcycle_stop(wcycle, ewcDOMDEC);
677 static void evaluate_energy(FILE *fplog, t_commrec *cr,
678 gmx_mtop_t *top_global,
679 em_state_t *ems, gmx_localtop_t *top,
680 t_inputrec *inputrec,
681 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
682 gmx_global_stat_t gstat,
683 gmx_vsite_t *vsite, gmx_constr_t constr,
685 t_graph *graph, t_mdatoms *mdatoms,
686 t_forcerec *fr, rvec mu_tot,
687 gmx_enerdata_t *enerd, tensor vir, tensor pres,
688 gmx_int64_t count, gmx_bool bFirst)
693 tensor force_vir, shake_vir, ekin;
694 real dvdl_constr, prescorr, enercorr, dvdlcorr;
697 /* Set the time to the initial time, the time does not change during EM */
698 t = inputrec->init_t;
701 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
703 /* This is the first state or an old state used before the last ns */
709 if (inputrec->nstlist > 0)
713 else if (inputrec->nstlist == -1)
715 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
718 gmx_sumi(1, &nabnsb, cr);
726 construct_vsites(vsite, ems->s.x, 1, NULL,
727 top->idef.iparams, top->idef.il,
728 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
731 if (DOMAINDECOMP(cr) && bNS)
733 /* Repartition the domain decomposition */
734 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
735 ems, top, mdatoms, fr, vsite, constr,
739 /* Calc force & energy on new trial position */
740 /* do_force always puts the charge groups in the box and shifts again
741 * We do not unshift, so molecules are always whole in congrad.c
743 do_force(fplog, cr, inputrec,
744 count, nrnb, wcycle, top, &top_global->groups,
745 ems->s.box, ems->s.x, &ems->s.hist,
746 ems->f, force_vir, mdatoms, enerd, fcd,
747 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
748 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
749 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
750 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
752 /* Clear the unused shake virial and pressure */
753 clear_mat(shake_vir);
756 /* Communicate stuff when parallel */
757 if (PAR(cr) && inputrec->eI != eiNM)
759 wallcycle_start(wcycle, ewcMoveE);
761 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
762 inputrec, NULL, NULL, NULL, 1, &terminate,
763 top_global, &ems->s, FALSE,
769 wallcycle_stop(wcycle, ewcMoveE);
772 /* Calculate long range corrections to pressure and energy */
773 calc_dispcorr(inputrec, fr, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
774 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
775 enerd->term[F_DISPCORR] = enercorr;
776 enerd->term[F_EPOT] += enercorr;
777 enerd->term[F_PRES] += prescorr;
778 enerd->term[F_DVDL] += dvdlcorr;
780 ems->epot = enerd->term[F_EPOT];
784 /* Project out the constraint components of the force */
785 wallcycle_start(wcycle, ewcCONSTR);
787 constrain(NULL, FALSE, FALSE, constr, &top->idef,
788 inputrec, NULL, cr, count, 0, 1.0, mdatoms,
789 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
790 ems->s.lambda[efptBONDED], &dvdl_constr,
791 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
792 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
793 m_add(force_vir, shake_vir, vir);
794 wallcycle_stop(wcycle, ewcCONSTR);
798 copy_mat(force_vir, vir);
802 enerd->term[F_PRES] =
803 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
805 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
807 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
809 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
813 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
815 em_state_t *s_min, em_state_t *s_b)
819 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
821 unsigned char *grpnrFREEZE;
825 fprintf(debug, "Doing reorder_partsum\n");
831 cgs_gl = dd_charge_groups_global(cr->dd);
832 index = cgs_gl->index;
834 /* Collect fm in a global vector fmg.
835 * This conflicts with the spirit of domain decomposition,
836 * but to fully optimize this a much more complicated algorithm is required.
838 snew(fmg, mtop->natoms);
840 ncg = s_min->s.ncg_gl;
841 cg_gl = s_min->s.cg_gl;
843 for (c = 0; c < ncg; c++)
848 for (a = a0; a < a1; a++)
850 copy_rvec(fm[i], fmg[a]);
854 gmx_sum(mtop->natoms*3, fmg[0], cr);
856 /* Now we will determine the part of the sum for the cgs in state s_b */
858 cg_gl = s_b->s.cg_gl;
862 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
863 for (c = 0; c < ncg; c++)
868 for (a = a0; a < a1; a++)
870 if (mdatoms->cFREEZE && grpnrFREEZE)
874 for (m = 0; m < DIM; m++)
876 if (!opts->nFreeze[gf][m])
878 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
890 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
892 em_state_t *s_min, em_state_t *s_b)
898 /* This is just the classical Polak-Ribiere calculation of beta;
899 * it looks a bit complicated since we take freeze groups into account,
900 * and might have to sum it in parallel runs.
903 if (!DOMAINDECOMP(cr) ||
904 (s_min->s.ddp_count == cr->dd->ddp_count &&
905 s_b->s.ddp_count == cr->dd->ddp_count))
911 /* This part of code can be incorrect with DD,
912 * since the atom ordering in s_b and s_min might differ.
914 for (i = 0; i < mdatoms->homenr; i++)
916 if (mdatoms->cFREEZE)
918 gf = mdatoms->cFREEZE[i];
920 for (m = 0; m < DIM; m++)
922 if (!opts->nFreeze[gf][m])
924 sum += (fb[i][m] - fm[i][m])*fb[i][m];
931 /* We need to reorder cgs while summing */
932 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
936 gmx_sumd(1, &sum, cr);
939 return sum/sqr(s_min->fnorm);
942 double do_cg(FILE *fplog, t_commrec *cr,
943 int nfile, const t_filenm fnm[],
944 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
945 int gmx_unused nstglobalcomm,
946 gmx_vsite_t *vsite, gmx_constr_t constr,
947 int gmx_unused stepout,
948 t_inputrec *inputrec,
949 gmx_mtop_t *top_global, t_fcdata *fcd,
950 t_state *state_global,
952 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
953 gmx_edsam_t gmx_unused ed,
955 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
956 gmx_membed_t gmx_unused membed,
957 real gmx_unused cpt_period, real gmx_unused max_hours,
958 const char gmx_unused *deviceOptions,
960 unsigned long gmx_unused Flags,
961 gmx_walltime_accounting_t walltime_accounting)
963 const char *CG = "Polak-Ribiere Conjugate Gradients";
965 em_state_t *s_min, *s_a, *s_b, *s_c;
967 gmx_enerdata_t *enerd;
969 gmx_global_stat_t gstat;
971 rvec *f_global, *p, *sf, *sfm;
972 double gpa, gpb, gpc, tmp, sum[2], minstep;
975 real a, b, c, beta = 0.0;
979 gmx_bool converged, foundlower;
981 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
983 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
985 int i, m, gf, step, nminstep;
990 s_min = init_em_state();
991 s_a = init_em_state();
992 s_b = init_em_state();
993 s_c = init_em_state();
995 /* Init em and store the local state in s_min */
996 init_em(fplog, CG, cr, inputrec,
997 state_global, top_global, s_min, &top, &f, &f_global,
998 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
999 nfile, fnm, &outf, &mdebin, imdport, Flags);
1001 /* Print to log file */
1002 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1004 /* Max number of steps */
1005 number_steps = inputrec->nsteps;
1009 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1013 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1016 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1017 /* do_force always puts the charge groups in the box and shifts again
1018 * We do not unshift, so molecules are always whole in congrad.c
1020 evaluate_energy(fplog, cr,
1021 top_global, s_min, top,
1022 inputrec, nrnb, wcycle, gstat,
1023 vsite, constr, fcd, graph, mdatoms, fr,
1024 mu_tot, enerd, vir, pres, -1, TRUE);
1029 /* Copy stuff to the energy bin for easy printing etc. */
1030 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1031 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1032 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1034 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1035 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1036 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1040 /* Estimate/guess the initial stepsize */
1041 stepsize = inputrec->em_stepsize/s_min->fnorm;
1045 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1046 s_min->fmax, s_min->a_fmax+1);
1047 fprintf(stderr, " F-Norm = %12.5e\n",
1048 s_min->fnorm/sqrt(state_global->natoms));
1049 fprintf(stderr, "\n");
1050 /* and copy to the log file too... */
1051 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1052 s_min->fmax, s_min->a_fmax+1);
1053 fprintf(fplog, " F-Norm = %12.5e\n",
1054 s_min->fnorm/sqrt(state_global->natoms));
1055 fprintf(fplog, "\n");
1057 /* Start the loop over CG steps.
1058 * Each successful step is counted, and we continue until
1059 * we either converge or reach the max number of steps.
1062 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1065 /* start taking steps in a new direction
1066 * First time we enter the routine, beta=0, and the direction is
1067 * simply the negative gradient.
1070 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1075 for (i = 0; i < mdatoms->homenr; i++)
1077 if (mdatoms->cFREEZE)
1079 gf = mdatoms->cFREEZE[i];
1081 for (m = 0; m < DIM; m++)
1083 if (!inputrec->opts.nFreeze[gf][m])
1085 p[i][m] = sf[i][m] + beta*p[i][m];
1086 gpa -= p[i][m]*sf[i][m];
1087 /* f is negative gradient, thus the sign */
1096 /* Sum the gradient along the line across CPUs */
1099 gmx_sumd(1, &gpa, cr);
1102 /* Calculate the norm of the search vector */
1103 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1105 /* Just in case stepsize reaches zero due to numerical precision... */
1108 stepsize = inputrec->em_stepsize/pnorm;
1112 * Double check the value of the derivative in the search direction.
1113 * If it is positive it must be due to the old information in the
1114 * CG formula, so just remove that and start over with beta=0.
1115 * This corresponds to a steepest descent step.
1120 step--; /* Don't count this step since we are restarting */
1121 continue; /* Go back to the beginning of the big for-loop */
1124 /* Calculate minimum allowed stepsize, before the average (norm)
1125 * relative change in coordinate is smaller than precision
1128 for (i = 0; i < mdatoms->homenr; i++)
1130 for (m = 0; m < DIM; m++)
1132 tmp = fabs(s_min->s.x[i][m]);
1141 /* Add up from all CPUs */
1144 gmx_sumd(1, &minstep, cr);
1147 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1149 if (stepsize < minstep)
1155 /* Write coordinates if necessary */
1156 do_x = do_per_step(step, inputrec->nstxout);
1157 do_f = do_per_step(step, inputrec->nstfout);
1159 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1160 top_global, inputrec, step,
1161 s_min, state_global, f_global);
1163 /* Take a step downhill.
1164 * In theory, we should minimize the function along this direction.
1165 * That is quite possible, but it turns out to take 5-10 function evaluations
1166 * for each line. However, we dont really need to find the exact minimum -
1167 * it is much better to start a new CG step in a modified direction as soon
1168 * as we are close to it. This will save a lot of energy evaluations.
1170 * In practice, we just try to take a single step.
1171 * If it worked (i.e. lowered the energy), we increase the stepsize but
1172 * the continue straight to the next CG step without trying to find any minimum.
1173 * If it didn't work (higher energy), there must be a minimum somewhere between
1174 * the old position and the new one.
1176 * Due to the finite numerical accuracy, it turns out that it is a good idea
1177 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1178 * This leads to lower final energies in the tests I've done. / Erik
1180 s_a->epot = s_min->epot;
1182 c = a + stepsize; /* reference position along line is zero */
1184 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1186 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1187 s_min, top, mdatoms, fr, vsite, constr,
1191 /* Take a trial step (new coords in s_c) */
1192 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1193 constr, top, nrnb, wcycle, -1);
1196 /* Calculate energy for the trial step */
1197 evaluate_energy(fplog, cr,
1198 top_global, s_c, top,
1199 inputrec, nrnb, wcycle, gstat,
1200 vsite, constr, fcd, graph, mdatoms, fr,
1201 mu_tot, enerd, vir, pres, -1, FALSE);
1203 /* Calc derivative along line */
1207 for (i = 0; i < mdatoms->homenr; i++)
1209 for (m = 0; m < DIM; m++)
1211 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1214 /* Sum the gradient along the line across CPUs */
1217 gmx_sumd(1, &gpc, cr);
1220 /* This is the max amount of increase in energy we tolerate */
1221 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1223 /* Accept the step if the energy is lower, or if it is not significantly higher
1224 * and the line derivative is still negative.
1226 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1229 /* Great, we found a better energy. Increase step for next iteration
1230 * if we are still going down, decrease it otherwise
1234 stepsize *= 1.618034; /* The golden section */
1238 stepsize *= 0.618034; /* 1/golden section */
1243 /* New energy is the same or higher. We will have to do some work
1244 * to find a smaller value in the interval. Take smaller step next time!
1247 stepsize *= 0.618034;
1253 /* OK, if we didn't find a lower value we will have to locate one now - there must
1254 * be one in the interval [a=0,c].
1255 * The same thing is valid here, though: Don't spend dozens of iterations to find
1256 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1257 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1259 * I also have a safeguard for potentially really patological functions so we never
1260 * take more than 20 steps before we give up ...
1262 * If we already found a lower value we just skip this step and continue to the update.
1270 /* Select a new trial point.
1271 * If the derivatives at points a & c have different sign we interpolate to zero,
1272 * otherwise just do a bisection.
1274 if (gpa < 0 && gpc > 0)
1276 b = a + gpa*(a-c)/(gpc-gpa);
1283 /* safeguard if interpolation close to machine accuracy causes errors:
1284 * never go outside the interval
1286 if (b <= a || b >= c)
1291 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1293 /* Reload the old state */
1294 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1295 s_min, top, mdatoms, fr, vsite, constr,
1299 /* Take a trial step to this new point - new coords in s_b */
1300 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1301 constr, top, nrnb, wcycle, -1);
1304 /* Calculate energy for the trial step */
1305 evaluate_energy(fplog, cr,
1306 top_global, s_b, top,
1307 inputrec, nrnb, wcycle, gstat,
1308 vsite, constr, fcd, graph, mdatoms, fr,
1309 mu_tot, enerd, vir, pres, -1, FALSE);
1311 /* p does not change within a step, but since the domain decomposition
1312 * might change, we have to use cg_p of s_b here.
1317 for (i = 0; i < mdatoms->homenr; i++)
1319 for (m = 0; m < DIM; m++)
1321 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1324 /* Sum the gradient along the line across CPUs */
1327 gmx_sumd(1, &gpb, cr);
1332 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1333 s_a->epot, s_b->epot, s_c->epot, gpb);
1336 epot_repl = s_b->epot;
1338 /* Keep one of the intervals based on the value of the derivative at the new point */
1341 /* Replace c endpoint with b */
1342 swap_em_state(s_b, s_c);
1348 /* Replace a endpoint with b */
1349 swap_em_state(s_b, s_a);
1355 * Stop search as soon as we find a value smaller than the endpoints.
1356 * Never run more than 20 steps, no matter what.
1360 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1363 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1366 /* OK. We couldn't find a significantly lower energy.
1367 * If beta==0 this was steepest descent, and then we give up.
1368 * If not, set beta=0 and restart with steepest descent before quitting.
1378 /* Reset memory before giving up */
1384 /* Select min energy state of A & C, put the best in B.
1386 if (s_c->epot < s_a->epot)
1390 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1391 s_c->epot, s_a->epot);
1393 swap_em_state(s_b, s_c);
1401 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1402 s_a->epot, s_c->epot);
1404 swap_em_state(s_b, s_a);
1414 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1417 swap_em_state(s_b, s_c);
1422 /* new search direction */
1423 /* beta = 0 means forget all memory and restart with steepest descents. */
1424 if (nstcg && ((step % nstcg) == 0))
1430 /* s_min->fnorm cannot be zero, because then we would have converged
1434 /* Polak-Ribiere update.
1435 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1437 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1439 /* Limit beta to prevent oscillations */
1440 if (fabs(beta) > 5.0)
1446 /* update positions */
1447 swap_em_state(s_min, s_b);
1450 /* Print it if necessary */
1455 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1456 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1457 s_min->fmax, s_min->a_fmax+1);
1459 /* Store the new (lower) energies */
1460 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1461 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1462 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1464 do_log = do_per_step(step, inputrec->nstlog);
1465 do_ene = do_per_step(step, inputrec->nstenergy);
1467 /* Prepare IMD energy record, if bIMD is TRUE. */
1468 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1472 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1474 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1475 do_log ? fplog : NULL, step, step, eprNORMAL,
1476 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1479 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1480 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1482 IMD_send_positions(inputrec->imd);
1485 /* Stop when the maximum force lies below tolerance.
1486 * If we have reached machine precision, converged is already set to true.
1488 converged = converged || (s_min->fmax < inputrec->em_tol);
1490 } /* End of the loop */
1492 /* IMD cleanup, if bIMD is TRUE. */
1493 IMD_finalize(inputrec->bIMD, inputrec->imd);
1497 step--; /* we never took that last step in this case */
1500 if (s_min->fmax > inputrec->em_tol)
1504 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1505 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1512 /* If we printed energy and/or logfile last step (which was the last step)
1513 * we don't have to do it again, but otherwise print the final values.
1517 /* Write final value to log since we didn't do anything the last step */
1518 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1520 if (!do_ene || !do_log)
1522 /* Write final energy file entries */
1523 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1524 !do_log ? fplog : NULL, step, step, eprNORMAL,
1525 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1529 /* Print some stuff... */
1532 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1536 * For accurate normal mode calculation it is imperative that we
1537 * store the last conformation into the full precision binary trajectory.
1539 * However, we should only do it if we did NOT already write this step
1540 * above (which we did if do_x or do_f was true).
1542 do_x = !do_per_step(step, inputrec->nstxout);
1543 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1545 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1546 top_global, inputrec, step,
1547 s_min, state_global, f_global);
1549 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1553 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1554 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1555 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1556 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1558 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1561 finish_em(cr, outf, walltime_accounting, wcycle);
1563 /* To print the actual number of steps we needed somewhere */
1564 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1567 } /* That's all folks */
1570 double do_lbfgs(FILE *fplog, t_commrec *cr,
1571 int nfile, const t_filenm fnm[],
1572 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1573 int gmx_unused nstglobalcomm,
1574 gmx_vsite_t *vsite, gmx_constr_t constr,
1575 int gmx_unused stepout,
1576 t_inputrec *inputrec,
1577 gmx_mtop_t *top_global, t_fcdata *fcd,
1580 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1581 gmx_edsam_t gmx_unused ed,
1583 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1584 gmx_membed_t gmx_unused membed,
1585 real gmx_unused cpt_period, real gmx_unused max_hours,
1586 const char gmx_unused *deviceOptions,
1588 unsigned long gmx_unused Flags,
1589 gmx_walltime_accounting_t walltime_accounting)
1591 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1593 gmx_localtop_t *top;
1594 gmx_enerdata_t *enerd;
1596 gmx_global_stat_t gstat;
1599 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1600 double stepsize, gpa, gpb, gpc, tmp, minstep;
1601 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1602 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1603 real a, b, c, maxdelta, delta;
1604 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1605 real dgdx, dgdg, sq, yr, beta;
1607 gmx_bool converged, first;
1610 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1612 int start, end, number_steps;
1614 int i, k, m, n, nfmax, gf, step;
1621 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1626 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).");
1629 n = 3*state->natoms;
1630 nmaxcorr = inputrec->nbfgscorr;
1632 /* Allocate memory */
1633 /* Use pointers to real so we dont have to loop over both atoms and
1634 * dimensions all the time...
1635 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1636 * that point to the same memory.
1649 snew(rho, nmaxcorr);
1650 snew(alpha, nmaxcorr);
1653 for (i = 0; i < nmaxcorr; i++)
1659 for (i = 0; i < nmaxcorr; i++)
1668 init_em(fplog, LBFGS, cr, inputrec,
1669 state, top_global, &ems, &top, &f, &f_global,
1670 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1671 nfile, fnm, &outf, &mdebin, imdport, Flags);
1672 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1673 * so we free some memory again.
1678 xx = (real *)state->x;
1682 end = mdatoms->homenr;
1684 /* Print to log file */
1685 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1687 do_log = do_ene = do_x = do_f = TRUE;
1689 /* Max number of steps */
1690 number_steps = inputrec->nsteps;
1692 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1694 for (i = start; i < end; i++)
1696 if (mdatoms->cFREEZE)
1698 gf = mdatoms->cFREEZE[i];
1700 for (m = 0; m < DIM; m++)
1702 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1707 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1711 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1716 construct_vsites(vsite, state->x, 1, NULL,
1717 top->idef.iparams, top->idef.il,
1718 fr->ePBC, fr->bMolPBC, cr, state->box);
1721 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1722 /* do_force always puts the charge groups in the box and shifts again
1723 * We do not unshift, so molecules are always whole
1728 evaluate_energy(fplog, cr,
1729 top_global, &ems, top,
1730 inputrec, nrnb, wcycle, gstat,
1731 vsite, constr, fcd, graph, mdatoms, fr,
1732 mu_tot, enerd, vir, pres, -1, TRUE);
1737 /* Copy stuff to the energy bin for easy printing etc. */
1738 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1739 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1740 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1742 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1743 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1744 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1748 /* This is the starting energy */
1749 Epot = enerd->term[F_EPOT];
1755 /* Set the initial step.
1756 * since it will be multiplied by the non-normalized search direction
1757 * vector (force vector the first time), we scale it by the
1758 * norm of the force.
1763 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1764 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1765 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1766 fprintf(stderr, "\n");
1767 /* and copy to the log file too... */
1768 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1769 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1770 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1771 fprintf(fplog, "\n");
1775 for (i = 0; i < n; i++)
1779 dx[point][i] = ff[i]; /* Initial search direction */
1787 stepsize = 1.0/fnorm;
1789 /* Start the loop over BFGS steps.
1790 * Each successful step is counted, and we continue until
1791 * we either converge or reach the max number of steps.
1796 /* Set the gradient from the force */
1798 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1801 /* Write coordinates if necessary */
1802 do_x = do_per_step(step, inputrec->nstxout);
1803 do_f = do_per_step(step, inputrec->nstfout);
1808 mdof_flags |= MDOF_X;
1813 mdof_flags |= MDOF_F;
1818 mdof_flags |= MDOF_IMD;
1821 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1822 top_global, step, (real)step, state, state, f, f);
1824 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1826 /* pointer to current direction - point=0 first time here */
1829 /* calculate line gradient */
1830 for (gpa = 0, i = 0; i < n; i++)
1835 /* Calculate minimum allowed stepsize, before the average (norm)
1836 * relative change in coordinate is smaller than precision
1838 for (minstep = 0, i = 0; i < n; i++)
1848 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1850 if (stepsize < minstep)
1856 /* Store old forces and coordinates */
1857 for (i = 0; i < n; i++)
1866 for (i = 0; i < n; i++)
1871 /* Take a step downhill.
1872 * In theory, we should minimize the function along this direction.
1873 * That is quite possible, but it turns out to take 5-10 function evaluations
1874 * for each line. However, we dont really need to find the exact minimum -
1875 * it is much better to start a new BFGS step in a modified direction as soon
1876 * as we are close to it. This will save a lot of energy evaluations.
1878 * In practice, we just try to take a single step.
1879 * If it worked (i.e. lowered the energy), we increase the stepsize but
1880 * the continue straight to the next BFGS step without trying to find any minimum.
1881 * If it didn't work (higher energy), there must be a minimum somewhere between
1882 * the old position and the new one.
1884 * Due to the finite numerical accuracy, it turns out that it is a good idea
1885 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1886 * This leads to lower final energies in the tests I've done. / Erik
1891 c = a + stepsize; /* reference position along line is zero */
1893 /* Check stepsize first. We do not allow displacements
1894 * larger than emstep.
1900 for (i = 0; i < n; i++)
1903 if (delta > maxdelta)
1908 if (maxdelta > inputrec->em_stepsize)
1913 while (maxdelta > inputrec->em_stepsize);
1915 /* Take a trial step */
1916 for (i = 0; i < n; i++)
1918 xc[i] = lastx[i] + c*s[i];
1922 /* Calculate energy for the trial step */
1923 ems.s.x = (rvec *)xc;
1925 evaluate_energy(fplog, cr,
1926 top_global, &ems, top,
1927 inputrec, nrnb, wcycle, gstat,
1928 vsite, constr, fcd, graph, mdatoms, fr,
1929 mu_tot, enerd, vir, pres, step, FALSE);
1932 /* Calc derivative along line */
1933 for (gpc = 0, i = 0; i < n; i++)
1935 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1937 /* Sum the gradient along the line across CPUs */
1940 gmx_sumd(1, &gpc, cr);
1943 /* This is the max amount of increase in energy we tolerate */
1944 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1946 /* Accept the step if the energy is lower, or if it is not significantly higher
1947 * and the line derivative is still negative.
1949 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1952 /* Great, we found a better energy. Increase step for next iteration
1953 * if we are still going down, decrease it otherwise
1957 stepsize *= 1.618034; /* The golden section */
1961 stepsize *= 0.618034; /* 1/golden section */
1966 /* New energy is the same or higher. We will have to do some work
1967 * to find a smaller value in the interval. Take smaller step next time!
1970 stepsize *= 0.618034;
1973 /* OK, if we didn't find a lower value we will have to locate one now - there must
1974 * be one in the interval [a=0,c].
1975 * The same thing is valid here, though: Don't spend dozens of iterations to find
1976 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1977 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1979 * I also have a safeguard for potentially really patological functions so we never
1980 * take more than 20 steps before we give up ...
1982 * If we already found a lower value we just skip this step and continue to the update.
1991 /* Select a new trial point.
1992 * If the derivatives at points a & c have different sign we interpolate to zero,
1993 * otherwise just do a bisection.
1996 if (gpa < 0 && gpc > 0)
1998 b = a + gpa*(a-c)/(gpc-gpa);
2005 /* safeguard if interpolation close to machine accuracy causes errors:
2006 * never go outside the interval
2008 if (b <= a || b >= c)
2013 /* Take a trial step */
2014 for (i = 0; i < n; i++)
2016 xb[i] = lastx[i] + b*s[i];
2020 /* Calculate energy for the trial step */
2021 ems.s.x = (rvec *)xb;
2023 evaluate_energy(fplog, cr,
2024 top_global, &ems, top,
2025 inputrec, nrnb, wcycle, gstat,
2026 vsite, constr, fcd, graph, mdatoms, fr,
2027 mu_tot, enerd, vir, pres, step, FALSE);
2032 for (gpb = 0, i = 0; i < n; i++)
2034 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2037 /* Sum the gradient along the line across CPUs */
2040 gmx_sumd(1, &gpb, cr);
2043 /* Keep one of the intervals based on the value of the derivative at the new point */
2046 /* Replace c endpoint with b */
2050 /* swap coord pointers b/c */
2060 /* Replace a endpoint with b */
2064 /* swap coord pointers a/b */
2074 * Stop search as soon as we find a value smaller than the endpoints,
2075 * or if the tolerance is below machine precision.
2076 * Never run more than 20 steps, no matter what.
2080 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2082 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2084 /* OK. We couldn't find a significantly lower energy.
2085 * If ncorr==0 this was steepest descent, and then we give up.
2086 * If not, reset memory to restart as steepest descent before quitting.
2098 /* Search in gradient direction */
2099 for (i = 0; i < n; i++)
2101 dx[point][i] = ff[i];
2103 /* Reset stepsize */
2104 stepsize = 1.0/fnorm;
2109 /* Select min energy state of A & C, put the best in xx/ff/Epot
2115 for (i = 0; i < n; i++)
2126 for (i = 0; i < n; i++)
2140 for (i = 0; i < n; i++)
2148 /* Update the memory information, and calculate a new
2149 * approximation of the inverse hessian
2152 /* Have new data in Epot, xx, ff */
2153 if (ncorr < nmaxcorr)
2158 for (i = 0; i < n; i++)
2160 dg[point][i] = lastf[i]-ff[i];
2161 dx[point][i] *= stepsize;
2166 for (i = 0; i < n; i++)
2168 dgdg += dg[point][i]*dg[point][i];
2169 dgdx += dg[point][i]*dx[point][i];
2174 rho[point] = 1.0/dgdx;
2177 if (point >= nmaxcorr)
2183 for (i = 0; i < n; i++)
2190 /* Recursive update. First go back over the memory points */
2191 for (k = 0; k < ncorr; k++)
2200 for (i = 0; i < n; i++)
2202 sq += dx[cp][i]*p[i];
2205 alpha[cp] = rho[cp]*sq;
2207 for (i = 0; i < n; i++)
2209 p[i] -= alpha[cp]*dg[cp][i];
2213 for (i = 0; i < n; i++)
2218 /* And then go forward again */
2219 for (k = 0; k < ncorr; k++)
2222 for (i = 0; i < n; i++)
2224 yr += p[i]*dg[cp][i];
2228 beta = alpha[cp]-beta;
2230 for (i = 0; i < n; i++)
2232 p[i] += beta*dx[cp][i];
2242 for (i = 0; i < n; i++)
2246 dx[point][i] = p[i];
2256 /* Test whether the convergence criterion is met */
2257 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2259 /* Print it if necessary */
2264 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2265 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2267 /* Store the new (lower) energies */
2268 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2269 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2270 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2271 do_log = do_per_step(step, inputrec->nstlog);
2272 do_ene = do_per_step(step, inputrec->nstenergy);
2275 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2277 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2278 do_log ? fplog : NULL, step, step, eprNORMAL,
2279 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2282 /* Send x and E to IMD client, if bIMD is TRUE. */
2283 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state->box, state->x, inputrec, 0, wcycle) && MASTER(cr))
2285 IMD_send_positions(inputrec->imd);
2288 /* Stop when the maximum force lies below tolerance.
2289 * If we have reached machine precision, converged is already set to true.
2292 converged = converged || (fmax < inputrec->em_tol);
2294 } /* End of the loop */
2296 /* IMD cleanup, if bIMD is TRUE. */
2297 IMD_finalize(inputrec->bIMD, inputrec->imd);
2301 step--; /* we never took that last step in this case */
2304 if (fmax > inputrec->em_tol)
2308 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2309 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2314 /* If we printed energy and/or logfile last step (which was the last step)
2315 * we don't have to do it again, but otherwise print the final values.
2317 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2319 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2321 if (!do_ene || !do_log) /* Write final energy file entries */
2323 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2324 !do_log ? fplog : NULL, step, step, eprNORMAL,
2325 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2328 /* Print some stuff... */
2331 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2335 * For accurate normal mode calculation it is imperative that we
2336 * store the last conformation into the full precision binary trajectory.
2338 * However, we should only do it if we did NOT already write this step
2339 * above (which we did if do_x or do_f was true).
2341 do_x = !do_per_step(step, inputrec->nstxout);
2342 do_f = !do_per_step(step, inputrec->nstfout);
2343 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2344 top_global, inputrec, step,
2349 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2350 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2351 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2352 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2354 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2357 finish_em(cr, outf, walltime_accounting, wcycle);
2359 /* To print the actual number of steps we needed somewhere */
2360 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2363 } /* That's all folks */
2366 double do_steep(FILE *fplog, t_commrec *cr,
2367 int nfile, const t_filenm fnm[],
2368 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2369 int gmx_unused nstglobalcomm,
2370 gmx_vsite_t *vsite, gmx_constr_t constr,
2371 int gmx_unused stepout,
2372 t_inputrec *inputrec,
2373 gmx_mtop_t *top_global, t_fcdata *fcd,
2374 t_state *state_global,
2376 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2377 gmx_edsam_t gmx_unused ed,
2379 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2380 gmx_membed_t gmx_unused membed,
2381 real gmx_unused cpt_period, real gmx_unused max_hours,
2382 const char gmx_unused *deviceOptions,
2384 unsigned long gmx_unused Flags,
2385 gmx_walltime_accounting_t walltime_accounting)
2387 const char *SD = "Steepest Descents";
2388 em_state_t *s_min, *s_try;
2390 gmx_localtop_t *top;
2391 gmx_enerdata_t *enerd;
2393 gmx_global_stat_t gstat;
2395 real stepsize, constepsize;
2399 gmx_bool bDone, bAbort, do_x, do_f;
2404 int steps_accepted = 0;
2408 s_min = init_em_state();
2409 s_try = init_em_state();
2411 /* Init em and store the local state in s_try */
2412 init_em(fplog, SD, cr, inputrec,
2413 state_global, top_global, s_try, &top, &f, &f_global,
2414 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2415 nfile, fnm, &outf, &mdebin, imdport, Flags);
2417 /* Print to log file */
2418 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2420 /* Set variables for stepsize (in nm). This is the largest
2421 * step that we are going to make in any direction.
2423 ustep = inputrec->em_stepsize;
2426 /* Max number of steps */
2427 nsteps = inputrec->nsteps;
2431 /* Print to the screen */
2432 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2436 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2439 /**** HERE STARTS THE LOOP ****
2440 * count is the counter for the number of steps
2441 * bDone will be TRUE when the minimization has converged
2442 * bAbort will be TRUE when nsteps steps have been performed or when
2443 * the stepsize becomes smaller than is reasonable for machine precision
2448 while (!bDone && !bAbort)
2450 bAbort = (nsteps >= 0) && (count == nsteps);
2452 /* set new coordinates, except for first step */
2455 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2456 s_min, stepsize, s_min->f, s_try,
2457 constr, top, nrnb, wcycle, count);
2460 evaluate_energy(fplog, cr,
2461 top_global, s_try, top,
2462 inputrec, nrnb, wcycle, gstat,
2463 vsite, constr, fcd, graph, mdatoms, fr,
2464 mu_tot, enerd, vir, pres, count, count == 0);
2468 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2473 s_min->epot = s_try->epot + 1;
2476 /* Print it if necessary */
2481 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2482 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2483 (s_try->epot < s_min->epot) ? '\n' : '\r');
2486 if (s_try->epot < s_min->epot)
2488 /* Store the new (lower) energies */
2489 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2490 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2491 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2493 /* Prepare IMD energy record, if bIMD is TRUE. */
2494 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2496 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2497 do_per_step(steps_accepted, inputrec->nstdisreout),
2498 do_per_step(steps_accepted, inputrec->nstorireout),
2499 fplog, count, count, eprNORMAL, TRUE,
2500 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2505 /* Now if the new energy is smaller than the previous...
2506 * or if this is the first step!
2507 * or if we did random steps!
2510 if ( (count == 0) || (s_try->epot < s_min->epot) )
2514 /* Test whether the convergence criterion is met... */
2515 bDone = (s_try->fmax < inputrec->em_tol);
2517 /* Copy the arrays for force, positions and energy */
2518 /* The 'Min' array always holds the coords and forces of the minimal
2520 swap_em_state(s_min, s_try);
2526 /* Write to trn, if necessary */
2527 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2528 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2529 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2530 top_global, inputrec, count,
2531 s_min, state_global, f_global);
2535 /* If energy is not smaller make the step smaller... */
2538 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2540 /* Reload the old state */
2541 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2542 s_min, top, mdatoms, fr, vsite, constr,
2547 /* Determine new step */
2548 stepsize = ustep/s_min->fmax;
2550 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2552 if (count == nsteps || ustep < 1e-12)
2554 if (count == nsteps || ustep < 1e-6)
2559 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2560 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2565 /* Send IMD energies and positions, if bIMD is TRUE. */
2566 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2568 IMD_send_positions(inputrec->imd);
2572 } /* End of the loop */
2574 /* IMD cleanup, if bIMD is TRUE. */
2575 IMD_finalize(inputrec->bIMD, inputrec->imd);
2577 /* Print some data... */
2580 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2582 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2583 top_global, inputrec, count,
2584 s_min, state_global, f_global);
2586 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2590 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2591 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2592 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2593 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2596 finish_em(cr, outf, walltime_accounting, wcycle);
2598 /* To print the actual number of steps we needed somewhere */
2599 inputrec->nsteps = count;
2601 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2604 } /* That's all folks */
2607 double do_nm(FILE *fplog, t_commrec *cr,
2608 int nfile, const t_filenm fnm[],
2609 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2610 int gmx_unused nstglobalcomm,
2611 gmx_vsite_t *vsite, gmx_constr_t constr,
2612 int gmx_unused stepout,
2613 t_inputrec *inputrec,
2614 gmx_mtop_t *top_global, t_fcdata *fcd,
2615 t_state *state_global,
2617 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2618 gmx_edsam_t gmx_unused ed,
2620 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2621 gmx_membed_t gmx_unused membed,
2622 real gmx_unused cpt_period, real gmx_unused max_hours,
2623 const char gmx_unused *deviceOptions,
2625 unsigned long gmx_unused Flags,
2626 gmx_walltime_accounting_t walltime_accounting)
2628 const char *NM = "Normal Mode Analysis";
2630 int natoms, atom, d;
2633 gmx_localtop_t *top;
2634 gmx_enerdata_t *enerd;
2636 gmx_global_stat_t gstat;
2638 real t, t0, lambda, lam0;
2643 gmx_bool bSparse; /* use sparse matrix storage format */
2645 gmx_sparsematrix_t * sparse_matrix = NULL;
2646 real * full_matrix = NULL;
2647 em_state_t * state_work;
2649 /* added with respect to mdrun */
2650 int i, j, k, row, col;
2651 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2657 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2660 state_work = init_em_state();
2662 /* Init em and store the local state in state_minimum */
2663 init_em(fplog, NM, cr, inputrec,
2664 state_global, top_global, state_work, &top,
2666 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2667 nfile, fnm, &outf, NULL, imdport, Flags);
2669 natoms = top_global->natoms;
2677 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2678 " which MIGHT not be accurate enough for normal mode analysis.\n"
2679 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2680 " are fairly modest even if you recompile in double precision.\n\n");
2684 /* Check if we can/should use sparse storage format.
2686 * Sparse format is only useful when the Hessian itself is sparse, which it
2687 * will be when we use a cutoff.
2688 * For small systems (n<1000) it is easier to always use full matrix format, though.
2690 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2692 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2695 else if (top_global->natoms < 1000)
2697 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2702 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2708 sz = DIM*top_global->natoms;
2710 fprintf(stderr, "Allocating Hessian memory...\n\n");
2714 sparse_matrix = gmx_sparsematrix_init(sz);
2715 sparse_matrix->compressed_symmetric = TRUE;
2719 snew(full_matrix, sz*sz);
2723 /* Initial values */
2724 t0 = inputrec->init_t;
2725 lam0 = inputrec->fepvals->init_lambda;
2733 /* Write start time and temperature */
2734 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2736 /* fudge nr of steps to nr of atoms */
2737 inputrec->nsteps = natoms*2;
2741 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2742 *(top_global->name), (int)inputrec->nsteps);
2745 nnodes = cr->nnodes;
2747 /* Make evaluate_energy do a single node force calculation */
2749 evaluate_energy(fplog, cr,
2750 top_global, state_work, top,
2751 inputrec, nrnb, wcycle, gstat,
2752 vsite, constr, fcd, graph, mdatoms, fr,
2753 mu_tot, enerd, vir, pres, -1, TRUE);
2754 cr->nnodes = nnodes;
2756 /* if forces are not small, warn user */
2757 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2759 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2760 if (state_work->fmax > 1.0e-3)
2762 md_print_info(cr, fplog,
2763 "The force is probably not small enough to "
2764 "ensure that you are at a minimum.\n"
2765 "Be aware that negative eigenvalues may occur\n"
2766 "when the resulting matrix is diagonalized.\n\n");
2769 /***********************************************************
2771 * Loop over all pairs in matrix
2773 * do_force called twice. Once with positive and
2774 * once with negative displacement
2776 ************************************************************/
2778 /* Steps are divided one by one over the nodes */
2779 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2782 for (d = 0; d < DIM; d++)
2784 x_min = state_work->s.x[atom][d];
2786 state_work->s.x[atom][d] = x_min - der_range;
2788 /* Make evaluate_energy do a single node force calculation */
2790 evaluate_energy(fplog, cr,
2791 top_global, state_work, top,
2792 inputrec, nrnb, wcycle, gstat,
2793 vsite, constr, fcd, graph, mdatoms, fr,
2794 mu_tot, enerd, vir, pres, atom*2, FALSE);
2796 for (i = 0; i < natoms; i++)
2798 copy_rvec(state_work->f[i], fneg[i]);
2801 state_work->s.x[atom][d] = x_min + der_range;
2803 evaluate_energy(fplog, cr,
2804 top_global, state_work, top,
2805 inputrec, nrnb, wcycle, gstat,
2806 vsite, constr, fcd, graph, mdatoms, fr,
2807 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2808 cr->nnodes = nnodes;
2810 /* x is restored to original */
2811 state_work->s.x[atom][d] = x_min;
2813 for (j = 0; j < natoms; j++)
2815 for (k = 0; (k < DIM); k++)
2818 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2826 #define mpi_type MPI_DOUBLE
2828 #define mpi_type MPI_FLOAT
2830 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2831 cr->mpi_comm_mygroup);
2836 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2842 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2843 cr->mpi_comm_mygroup, &stat);
2848 row = (atom + node)*DIM + d;
2850 for (j = 0; j < natoms; j++)
2852 for (k = 0; k < DIM; k++)
2858 if (col >= row && dfdx[j][k] != 0.0)
2860 gmx_sparsematrix_increment_value(sparse_matrix,
2861 row, col, dfdx[j][k]);
2866 full_matrix[row*sz+col] = dfdx[j][k];
2873 if (bVerbose && fplog)
2878 /* write progress */
2879 if (MASTER(cr) && bVerbose)
2881 fprintf(stderr, "\rFinished step %d out of %d",
2882 min(atom+nnodes, natoms), natoms);
2889 fprintf(stderr, "\n\nWriting Hessian...\n");
2890 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2893 finish_em(cr, outf, walltime_accounting, wcycle);
2895 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);