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.
59 #include "md_support.h"
64 #include "mtop_util.h"
67 #include "gmx_omp_nthreads.h"
68 #include "md_logging.h"
70 #include "gromacs/fileio/confio.h"
71 #include "gromacs/fileio/mtxio.h"
72 #include "gromacs/fileio/trajectory_writing.h"
73 #include "gromacs/imd/imd.h"
74 #include "gromacs/legacyheaders/types/commrec.h"
75 #include "gromacs/linearalgebra/sparsematrix.h"
76 #include "gromacs/math/vec.h"
77 #include "gromacs/pbcutil/mshift.h"
78 #include "gromacs/timing/wallcycle.h"
79 #include "gromacs/timing/walltime_accounting.h"
80 #include "gromacs/utility/cstringutil.h"
81 #include "gromacs/utility/fatalerror.h"
82 #include "gromacs/utility/smalloc.h"
93 static em_state_t *init_em_state()
99 /* does this need to be here? Should the array be declared differently (staticaly)in the state definition? */
100 snew(ems->s.lambda, efptNR);
105 static void print_em_start(FILE *fplog,
107 gmx_walltime_accounting_t walltime_accounting,
108 gmx_wallcycle_t wcycle,
111 walltime_accounting_start(walltime_accounting);
112 wallcycle_start(wcycle, ewcRUN);
113 print_start(fplog, cr, walltime_accounting, name);
115 static void em_time_end(gmx_walltime_accounting_t walltime_accounting,
116 gmx_wallcycle_t wcycle)
118 wallcycle_stop(wcycle, ewcRUN);
120 walltime_accounting_end(walltime_accounting);
123 static void sp_header(FILE *out, const char *minimizer, real ftol, int nsteps)
126 fprintf(out, "%s:\n", minimizer);
127 fprintf(out, " Tolerance (Fmax) = %12.5e\n", ftol);
128 fprintf(out, " Number of steps = %12d\n", nsteps);
131 static void warn_step(FILE *fp, real ftol, gmx_bool bLastStep, gmx_bool bConstrain)
137 "\nEnergy minimization reached the maximum number "
138 "of steps before the forces reached the requested "
139 "precision Fmax < %g.\n", ftol);
144 "\nEnergy minimization has stopped, but the forces have "
145 "not converged to the requested precision Fmax < %g (which "
146 "may not be possible for your system). It stopped "
147 "because the algorithm tried to make a new step whose size "
148 "was too small, or there was no change in the energy since "
149 "last step. Either way, we regard the minimization as "
150 "converged to within the available machine precision, "
151 "given your starting configuration and EM parameters.\n%s%s",
153 sizeof(real) < sizeof(double) ?
154 "\nDouble precision normally gives you higher accuracy, but "
155 "this is often not needed for preparing to run molecular "
159 "You might need to increase your constraint accuracy, or turn\n"
160 "off constraints altogether (set constraints = none in mdp file)\n" :
163 fputs(wrap_lines(buffer, 78, 0, FALSE), fp);
168 static void print_converged(FILE *fp, const char *alg, real ftol,
169 gmx_int64_t count, gmx_bool bDone, gmx_int64_t nsteps,
170 real epot, real fmax, int nfmax, real fnorm)
172 char buf[STEPSTRSIZE];
176 fprintf(fp, "\n%s converged to Fmax < %g in %s steps\n",
177 alg, ftol, gmx_step_str(count, buf));
179 else if (count < nsteps)
181 fprintf(fp, "\n%s converged to machine precision in %s steps,\n"
182 "but did not reach the requested Fmax < %g.\n",
183 alg, gmx_step_str(count, buf), ftol);
187 fprintf(fp, "\n%s did not converge to Fmax < %g in %s steps.\n",
188 alg, ftol, gmx_step_str(count, buf));
192 fprintf(fp, "Potential Energy = %21.14e\n", epot);
193 fprintf(fp, "Maximum force = %21.14e on atom %d\n", fmax, nfmax+1);
194 fprintf(fp, "Norm of force = %21.14e\n", fnorm);
196 fprintf(fp, "Potential Energy = %14.7e\n", epot);
197 fprintf(fp, "Maximum force = %14.7e on atom %d\n", fmax, nfmax+1);
198 fprintf(fp, "Norm of force = %14.7e\n", fnorm);
202 static void get_f_norm_max(t_commrec *cr,
203 t_grpopts *opts, t_mdatoms *mdatoms, rvec *f,
204 real *fnorm, real *fmax, int *a_fmax)
207 real fmax2, fmax2_0, fam;
208 int la_max, a_max, start, end, i, m, gf;
210 /* This routine finds the largest force and returns it.
211 * On parallel machines the global max is taken.
218 end = mdatoms->homenr;
219 if (mdatoms->cFREEZE)
221 for (i = start; i < end; i++)
223 gf = mdatoms->cFREEZE[i];
225 for (m = 0; m < DIM; m++)
227 if (!opts->nFreeze[gf][m])
242 for (i = start; i < end; i++)
254 if (la_max >= 0 && DOMAINDECOMP(cr))
256 a_max = cr->dd->gatindex[la_max];
264 snew(sum, 2*cr->nnodes+1);
265 sum[2*cr->nodeid] = fmax2;
266 sum[2*cr->nodeid+1] = a_max;
267 sum[2*cr->nnodes] = fnorm2;
268 gmx_sumd(2*cr->nnodes+1, sum, cr);
269 fnorm2 = sum[2*cr->nnodes];
270 /* Determine the global maximum */
271 for (i = 0; i < cr->nnodes; i++)
273 if (sum[2*i] > fmax2)
276 a_max = (int)(sum[2*i+1] + 0.5);
284 *fnorm = sqrt(fnorm2);
296 static void get_state_f_norm_max(t_commrec *cr,
297 t_grpopts *opts, t_mdatoms *mdatoms,
300 get_f_norm_max(cr, opts, mdatoms, ems->f, &ems->fnorm, &ems->fmax, &ems->a_fmax);
303 void init_em(FILE *fplog, const char *title,
304 t_commrec *cr, t_inputrec *ir,
305 t_state *state_global, gmx_mtop_t *top_global,
306 em_state_t *ems, gmx_localtop_t **top,
307 rvec **f, rvec **f_global,
308 t_nrnb *nrnb, rvec mu_tot,
309 t_forcerec *fr, gmx_enerdata_t **enerd,
310 t_graph **graph, t_mdatoms *mdatoms, gmx_global_stat_t *gstat,
311 gmx_vsite_t *vsite, gmx_constr_t constr,
312 int nfile, const t_filenm fnm[],
313 gmx_mdoutf_t *outf, t_mdebin **mdebin,
314 int imdport, unsigned long gmx_unused Flags)
321 fprintf(fplog, "Initiating %s\n", title);
324 state_global->ngtc = 0;
326 /* Initialize lambda variables */
327 initialize_lambdas(fplog, ir, &(state_global->fep_state), state_global->lambda, NULL);
331 /* Interactive molecular dynamics */
332 init_IMD(ir, cr, top_global, fplog, 1, state_global->x,
333 nfile, fnm, NULL, imdport, Flags);
335 if (DOMAINDECOMP(cr))
337 *top = dd_init_local_top(top_global);
339 dd_init_local_state(cr->dd, state_global, &ems->s);
343 /* Distribute the charge groups over the nodes from the master node */
344 dd_partition_system(fplog, ir->init_step, cr, TRUE, 1,
345 state_global, top_global, ir,
346 &ems->s, &ems->f, mdatoms, *top,
347 fr, vsite, NULL, constr,
349 dd_store_state(cr->dd, &ems->s);
353 snew(*f_global, top_global->natoms);
363 snew(*f, top_global->natoms);
365 /* Just copy the state */
366 ems->s = *state_global;
367 snew(ems->s.x, ems->s.nalloc);
368 snew(ems->f, ems->s.nalloc);
369 for (i = 0; i < state_global->natoms; i++)
371 copy_rvec(state_global->x[i], ems->s.x[i]);
373 copy_mat(state_global->box, ems->s.box);
375 *top = gmx_mtop_generate_local_top(top_global, ir);
378 forcerec_set_excl_load(fr, *top);
380 setup_bonded_threading(fr, &(*top)->idef);
382 if (ir->ePBC != epbcNONE && !fr->bMolPBC)
384 *graph = mk_graph(fplog, &((*top)->idef), 0, top_global->natoms, FALSE, FALSE);
391 atoms2md(top_global, ir, 0, NULL, top_global->natoms, mdatoms);
392 update_mdatoms(mdatoms, state_global->lambda[efptFEP]);
396 set_vsite_top(vsite, *top, mdatoms, cr);
402 if (ir->eConstrAlg == econtSHAKE &&
403 gmx_mtop_ftype_count(top_global, F_CONSTR) > 0)
405 gmx_fatal(FARGS, "Can not do energy minimization with %s, use %s\n",
406 econstr_names[econtSHAKE], econstr_names[econtLINCS]);
409 if (!DOMAINDECOMP(cr))
411 set_constraints(constr, *top, ir, mdatoms, cr);
414 if (!ir->bContinuation)
416 /* Constrain the starting coordinates */
418 constrain(PAR(cr) ? NULL : fplog, TRUE, TRUE, constr, &(*top)->idef,
419 ir, NULL, cr, -1, 0, mdatoms,
420 ems->s.x, ems->s.x, NULL, fr->bMolPBC, ems->s.box,
421 ems->s.lambda[efptFEP], &dvdl_constr,
422 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
428 *gstat = global_stat_init(ir);
431 *outf = init_mdoutf(fplog, nfile, fnm, 0, cr, ir, top_global, NULL);
434 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
439 /* Init bin for energy stuff */
440 *mdebin = init_mdebin(mdoutf_get_fp_ene(*outf), top_global, ir, NULL);
444 calc_shifts(ems->s.box, fr->shift_vec);
447 static void finish_em(t_commrec *cr, gmx_mdoutf_t outf,
448 gmx_walltime_accounting_t walltime_accounting,
449 gmx_wallcycle_t wcycle)
451 if (!(cr->duty & DUTY_PME))
453 /* Tell the PME only node to finish */
454 gmx_pme_send_finish(cr);
459 em_time_end(walltime_accounting, wcycle);
462 static void swap_em_state(em_state_t *ems1, em_state_t *ems2)
471 static void copy_em_coords(em_state_t *ems, t_state *state)
475 for (i = 0; (i < state->natoms); i++)
477 copy_rvec(ems->s.x[i], state->x[i]);
481 static void write_em_traj(FILE *fplog, t_commrec *cr,
483 gmx_bool bX, gmx_bool bF, const char *confout,
484 gmx_mtop_t *top_global,
485 t_inputrec *ir, gmx_int64_t step,
487 t_state *state_global, rvec *f_global)
490 gmx_bool bIMDout = FALSE;
493 /* Shall we do IMD output? */
496 bIMDout = do_per_step(step, IMD_get_step(ir->imd->setup));
499 if ((bX || bF || bIMDout || confout != NULL) && !DOMAINDECOMP(cr))
501 copy_em_coords(state, state_global);
508 mdof_flags |= MDOF_X;
512 mdof_flags |= MDOF_F;
515 /* If we want IMD output, set appropriate MDOF flag */
518 mdof_flags |= MDOF_IMD;
521 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
522 top_global, step, (double)step,
523 &state->s, state_global, state->f, f_global);
525 if (confout != NULL && MASTER(cr))
527 if (ir->ePBC != epbcNONE && !ir->bPeriodicMols && DOMAINDECOMP(cr))
529 /* Make molecules whole only for confout writing */
530 do_pbc_mtop(fplog, ir->ePBC, state_global->box, top_global,
534 write_sto_conf_mtop(confout,
535 *top_global->name, top_global,
536 state_global->x, NULL, ir->ePBC, state_global->box);
540 static void do_em_step(t_commrec *cr, t_inputrec *ir, t_mdatoms *md,
542 em_state_t *ems1, real a, rvec *f, em_state_t *ems2,
543 gmx_constr_t constr, gmx_localtop_t *top,
544 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
557 if (DOMAINDECOMP(cr) && s1->ddp_count != cr->dd->ddp_count)
559 gmx_incons("state mismatch in do_em_step");
562 s2->flags = s1->flags;
564 if (s2->nalloc != s1->nalloc)
566 s2->nalloc = s1->nalloc;
567 srenew(s2->x, s1->nalloc);
568 srenew(ems2->f, s1->nalloc);
569 if (s2->flags & (1<<estCGP))
571 srenew(s2->cg_p, s1->nalloc);
575 s2->natoms = s1->natoms;
576 copy_mat(s1->box, s2->box);
577 /* Copy free energy state */
578 for (i = 0; i < efptNR; i++)
580 s2->lambda[i] = s1->lambda[i];
582 copy_mat(s1->box, s2->box);
590 #pragma omp parallel num_threads(gmx_omp_nthreads_get(emntUpdate))
595 #pragma omp for schedule(static) nowait
596 for (i = start; i < end; i++)
602 for (m = 0; m < DIM; m++)
604 if (ir->opts.nFreeze[gf][m])
610 x2[i][m] = x1[i][m] + a*f[i][m];
615 if (s2->flags & (1<<estCGP))
617 /* Copy the CG p vector */
620 #pragma omp for schedule(static) nowait
621 for (i = start; i < end; i++)
623 copy_rvec(x1[i], x2[i]);
627 if (DOMAINDECOMP(cr))
629 s2->ddp_count = s1->ddp_count;
630 if (s2->cg_gl_nalloc < s1->cg_gl_nalloc)
633 s2->cg_gl_nalloc = s1->cg_gl_nalloc;
634 srenew(s2->cg_gl, s2->cg_gl_nalloc);
637 s2->ncg_gl = s1->ncg_gl;
638 #pragma omp for schedule(static) nowait
639 for (i = 0; i < s2->ncg_gl; i++)
641 s2->cg_gl[i] = s1->cg_gl[i];
643 s2->ddp_count_cg_gl = s1->ddp_count_cg_gl;
649 wallcycle_start(wcycle, ewcCONSTR);
651 constrain(NULL, TRUE, TRUE, constr, &top->idef,
652 ir, NULL, cr, count, 0, md,
653 s1->x, s2->x, NULL, bMolPBC, s2->box,
654 s2->lambda[efptBONDED], &dvdl_constr,
655 NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
656 wallcycle_stop(wcycle, ewcCONSTR);
660 static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
661 gmx_mtop_t *top_global, t_inputrec *ir,
662 em_state_t *ems, gmx_localtop_t *top,
663 t_mdatoms *mdatoms, t_forcerec *fr,
664 gmx_vsite_t *vsite, gmx_constr_t constr,
665 t_nrnb *nrnb, gmx_wallcycle_t wcycle)
667 /* Repartition the domain decomposition */
668 wallcycle_start(wcycle, ewcDOMDEC);
669 dd_partition_system(fplog, step, cr, FALSE, 1,
670 NULL, top_global, ir,
672 mdatoms, top, fr, vsite, NULL, constr,
673 nrnb, wcycle, FALSE);
674 dd_store_state(cr->dd, &ems->s);
675 wallcycle_stop(wcycle, ewcDOMDEC);
678 static void evaluate_energy(FILE *fplog, t_commrec *cr,
679 gmx_mtop_t *top_global,
680 em_state_t *ems, gmx_localtop_t *top,
681 t_inputrec *inputrec,
682 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
683 gmx_global_stat_t gstat,
684 gmx_vsite_t *vsite, gmx_constr_t constr,
686 t_graph *graph, t_mdatoms *mdatoms,
687 t_forcerec *fr, rvec mu_tot,
688 gmx_enerdata_t *enerd, tensor vir, tensor pres,
689 gmx_int64_t count, gmx_bool bFirst)
694 tensor force_vir, shake_vir, ekin;
695 real dvdl_constr, prescorr, enercorr, dvdlcorr;
698 /* Set the time to the initial time, the time does not change during EM */
699 t = inputrec->init_t;
702 (DOMAINDECOMP(cr) && ems->s.ddp_count < cr->dd->ddp_count))
704 /* This is the first state or an old state used before the last ns */
710 if (inputrec->nstlist > 0)
714 else if (inputrec->nstlist == -1)
716 nabnsb = natoms_beyond_ns_buffer(inputrec, fr, &top->cgs, NULL, ems->s.x);
719 gmx_sumi(1, &nabnsb, cr);
727 construct_vsites(vsite, ems->s.x, 1, NULL,
728 top->idef.iparams, top->idef.il,
729 fr->ePBC, fr->bMolPBC, cr, ems->s.box);
732 if (DOMAINDECOMP(cr) && bNS)
734 /* Repartition the domain decomposition */
735 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
736 ems, top, mdatoms, fr, vsite, constr,
740 /* Calc force & energy on new trial position */
741 /* do_force always puts the charge groups in the box and shifts again
742 * We do not unshift, so molecules are always whole in congrad.c
744 do_force(fplog, cr, inputrec,
745 count, nrnb, wcycle, top, &top_global->groups,
746 ems->s.box, ems->s.x, &ems->s.hist,
747 ems->f, force_vir, mdatoms, enerd, fcd,
748 ems->s.lambda, graph, fr, vsite, mu_tot, t, NULL, NULL, TRUE,
749 GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES |
750 GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY |
751 (bNS ? GMX_FORCE_NS | GMX_FORCE_DO_LR : 0));
753 /* Clear the unused shake virial and pressure */
754 clear_mat(shake_vir);
757 /* Communicate stuff when parallel */
758 if (PAR(cr) && inputrec->eI != eiNM)
760 wallcycle_start(wcycle, ewcMoveE);
762 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
763 inputrec, NULL, NULL, NULL, 1, &terminate,
764 top_global, &ems->s, FALSE,
770 wallcycle_stop(wcycle, ewcMoveE);
773 /* Calculate long range corrections to pressure and energy */
774 calc_dispcorr(fplog, inputrec, fr, count, top_global->natoms, ems->s.box, ems->s.lambda[efptVDW],
775 pres, force_vir, &prescorr, &enercorr, &dvdlcorr);
776 enerd->term[F_DISPCORR] = enercorr;
777 enerd->term[F_EPOT] += enercorr;
778 enerd->term[F_PRES] += prescorr;
779 enerd->term[F_DVDL] += dvdlcorr;
781 ems->epot = enerd->term[F_EPOT];
785 /* Project out the constraint components of the force */
786 wallcycle_start(wcycle, ewcCONSTR);
788 constrain(NULL, FALSE, FALSE, constr, &top->idef,
789 inputrec, NULL, cr, count, 0, mdatoms,
790 ems->s.x, ems->f, ems->f, fr->bMolPBC, ems->s.box,
791 ems->s.lambda[efptBONDED], &dvdl_constr,
792 NULL, &shake_vir, nrnb, econqForceDispl, FALSE, 0, 0);
793 if (fr->bSepDVDL && fplog)
795 gmx_print_sepdvdl(fplog, "Constraints", t, dvdl_constr);
797 enerd->term[F_DVDL_CONSTR] += dvdl_constr;
798 m_add(force_vir, shake_vir, vir);
799 wallcycle_stop(wcycle, ewcCONSTR);
803 copy_mat(force_vir, vir);
807 enerd->term[F_PRES] =
808 calc_pres(fr->ePBC, inputrec->nwall, ems->s.box, ekin, vir, pres);
810 sum_dhdl(enerd, ems->s.lambda, inputrec->fepvals);
812 if (EI_ENERGY_MINIMIZATION(inputrec->eI))
814 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, ems);
818 static double reorder_partsum(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
820 em_state_t *s_min, em_state_t *s_b)
824 int ncg, *cg_gl, *index, c, cg, i, a0, a1, a, gf, m;
826 unsigned char *grpnrFREEZE;
830 fprintf(debug, "Doing reorder_partsum\n");
836 cgs_gl = dd_charge_groups_global(cr->dd);
837 index = cgs_gl->index;
839 /* Collect fm in a global vector fmg.
840 * This conflicts with the spirit of domain decomposition,
841 * but to fully optimize this a much more complicated algorithm is required.
843 snew(fmg, mtop->natoms);
845 ncg = s_min->s.ncg_gl;
846 cg_gl = s_min->s.cg_gl;
848 for (c = 0; c < ncg; c++)
853 for (a = a0; a < a1; a++)
855 copy_rvec(fm[i], fmg[a]);
859 gmx_sum(mtop->natoms*3, fmg[0], cr);
861 /* Now we will determine the part of the sum for the cgs in state s_b */
863 cg_gl = s_b->s.cg_gl;
867 grpnrFREEZE = mtop->groups.grpnr[egcFREEZE];
868 for (c = 0; c < ncg; c++)
873 for (a = a0; a < a1; a++)
875 if (mdatoms->cFREEZE && grpnrFREEZE)
879 for (m = 0; m < DIM; m++)
881 if (!opts->nFreeze[gf][m])
883 partsum += (fb[i][m] - fmg[a][m])*fb[i][m];
895 static real pr_beta(t_commrec *cr, t_grpopts *opts, t_mdatoms *mdatoms,
897 em_state_t *s_min, em_state_t *s_b)
903 /* This is just the classical Polak-Ribiere calculation of beta;
904 * it looks a bit complicated since we take freeze groups into account,
905 * and might have to sum it in parallel runs.
908 if (!DOMAINDECOMP(cr) ||
909 (s_min->s.ddp_count == cr->dd->ddp_count &&
910 s_b->s.ddp_count == cr->dd->ddp_count))
916 /* This part of code can be incorrect with DD,
917 * since the atom ordering in s_b and s_min might differ.
919 for (i = 0; i < mdatoms->homenr; i++)
921 if (mdatoms->cFREEZE)
923 gf = mdatoms->cFREEZE[i];
925 for (m = 0; m < DIM; m++)
927 if (!opts->nFreeze[gf][m])
929 sum += (fb[i][m] - fm[i][m])*fb[i][m];
936 /* We need to reorder cgs while summing */
937 sum = reorder_partsum(cr, opts, mdatoms, mtop, s_min, s_b);
941 gmx_sumd(1, &sum, cr);
944 return sum/sqr(s_min->fnorm);
947 double do_cg(FILE *fplog, t_commrec *cr,
948 int nfile, const t_filenm fnm[],
949 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
950 int gmx_unused nstglobalcomm,
951 gmx_vsite_t *vsite, gmx_constr_t constr,
952 int gmx_unused stepout,
953 t_inputrec *inputrec,
954 gmx_mtop_t *top_global, t_fcdata *fcd,
955 t_state *state_global,
957 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
958 gmx_edsam_t gmx_unused ed,
960 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
961 gmx_membed_t gmx_unused membed,
962 real gmx_unused cpt_period, real gmx_unused max_hours,
963 const char gmx_unused *deviceOptions,
965 unsigned long gmx_unused Flags,
966 gmx_walltime_accounting_t walltime_accounting)
968 const char *CG = "Polak-Ribiere Conjugate Gradients";
970 em_state_t *s_min, *s_a, *s_b, *s_c;
972 gmx_enerdata_t *enerd;
974 gmx_global_stat_t gstat;
976 rvec *f_global, *p, *sf, *sfm;
977 double gpa, gpb, gpc, tmp, sum[2], minstep;
980 real a, b, c, beta = 0.0;
984 gmx_bool converged, foundlower;
986 gmx_bool do_log = FALSE, do_ene = FALSE, do_x, do_f;
988 int number_steps, neval = 0, nstcg = inputrec->nstcgsteep;
990 int i, m, gf, step, nminstep;
995 s_min = init_em_state();
996 s_a = init_em_state();
997 s_b = init_em_state();
998 s_c = init_em_state();
1000 /* Init em and store the local state in s_min */
1001 init_em(fplog, CG, cr, inputrec,
1002 state_global, top_global, s_min, &top, &f, &f_global,
1003 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1004 nfile, fnm, &outf, &mdebin, imdport, Flags);
1006 /* Print to log file */
1007 print_em_start(fplog, cr, walltime_accounting, wcycle, CG);
1009 /* Max number of steps */
1010 number_steps = inputrec->nsteps;
1014 sp_header(stderr, CG, inputrec->em_tol, number_steps);
1018 sp_header(fplog, CG, inputrec->em_tol, number_steps);
1021 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1022 /* do_force always puts the charge groups in the box and shifts again
1023 * We do not unshift, so molecules are always whole in congrad.c
1025 evaluate_energy(fplog, cr,
1026 top_global, s_min, top,
1027 inputrec, nrnb, wcycle, gstat,
1028 vsite, constr, fcd, graph, mdatoms, fr,
1029 mu_tot, enerd, vir, pres, -1, TRUE);
1034 /* Copy stuff to the energy bin for easy printing etc. */
1035 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1036 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1037 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1039 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1040 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1041 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1045 /* Estimate/guess the initial stepsize */
1046 stepsize = inputrec->em_stepsize/s_min->fnorm;
1050 fprintf(stderr, " F-max = %12.5e on atom %d\n",
1051 s_min->fmax, s_min->a_fmax+1);
1052 fprintf(stderr, " F-Norm = %12.5e\n",
1053 s_min->fnorm/sqrt(state_global->natoms));
1054 fprintf(stderr, "\n");
1055 /* and copy to the log file too... */
1056 fprintf(fplog, " F-max = %12.5e on atom %d\n",
1057 s_min->fmax, s_min->a_fmax+1);
1058 fprintf(fplog, " F-Norm = %12.5e\n",
1059 s_min->fnorm/sqrt(state_global->natoms));
1060 fprintf(fplog, "\n");
1062 /* Start the loop over CG steps.
1063 * Each successful step is counted, and we continue until
1064 * we either converge or reach the max number of steps.
1067 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1070 /* start taking steps in a new direction
1071 * First time we enter the routine, beta=0, and the direction is
1072 * simply the negative gradient.
1075 /* Calculate the new direction in p, and the gradient in this direction, gpa */
1080 for (i = 0; i < mdatoms->homenr; i++)
1082 if (mdatoms->cFREEZE)
1084 gf = mdatoms->cFREEZE[i];
1086 for (m = 0; m < DIM; m++)
1088 if (!inputrec->opts.nFreeze[gf][m])
1090 p[i][m] = sf[i][m] + beta*p[i][m];
1091 gpa -= p[i][m]*sf[i][m];
1092 /* f is negative gradient, thus the sign */
1101 /* Sum the gradient along the line across CPUs */
1104 gmx_sumd(1, &gpa, cr);
1107 /* Calculate the norm of the search vector */
1108 get_f_norm_max(cr, &(inputrec->opts), mdatoms, p, &pnorm, NULL, NULL);
1110 /* Just in case stepsize reaches zero due to numerical precision... */
1113 stepsize = inputrec->em_stepsize/pnorm;
1117 * Double check the value of the derivative in the search direction.
1118 * If it is positive it must be due to the old information in the
1119 * CG formula, so just remove that and start over with beta=0.
1120 * This corresponds to a steepest descent step.
1125 step--; /* Don't count this step since we are restarting */
1126 continue; /* Go back to the beginning of the big for-loop */
1129 /* Calculate minimum allowed stepsize, before the average (norm)
1130 * relative change in coordinate is smaller than precision
1133 for (i = 0; i < mdatoms->homenr; i++)
1135 for (m = 0; m < DIM; m++)
1137 tmp = fabs(s_min->s.x[i][m]);
1146 /* Add up from all CPUs */
1149 gmx_sumd(1, &minstep, cr);
1152 minstep = GMX_REAL_EPS/sqrt(minstep/(3*state_global->natoms));
1154 if (stepsize < minstep)
1160 /* Write coordinates if necessary */
1161 do_x = do_per_step(step, inputrec->nstxout);
1162 do_f = do_per_step(step, inputrec->nstfout);
1164 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
1165 top_global, inputrec, step,
1166 s_min, state_global, f_global);
1168 /* Take a step downhill.
1169 * In theory, we should minimize the function along this direction.
1170 * That is quite possible, but it turns out to take 5-10 function evaluations
1171 * for each line. However, we dont really need to find the exact minimum -
1172 * it is much better to start a new CG step in a modified direction as soon
1173 * as we are close to it. This will save a lot of energy evaluations.
1175 * In practice, we just try to take a single step.
1176 * If it worked (i.e. lowered the energy), we increase the stepsize but
1177 * the continue straight to the next CG step without trying to find any minimum.
1178 * If it didn't work (higher energy), there must be a minimum somewhere between
1179 * the old position and the new one.
1181 * Due to the finite numerical accuracy, it turns out that it is a good idea
1182 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1183 * This leads to lower final energies in the tests I've done. / Erik
1185 s_a->epot = s_min->epot;
1187 c = a + stepsize; /* reference position along line is zero */
1189 if (DOMAINDECOMP(cr) && s_min->s.ddp_count < cr->dd->ddp_count)
1191 em_dd_partition_system(fplog, step, cr, top_global, inputrec,
1192 s_min, top, mdatoms, fr, vsite, constr,
1196 /* Take a trial step (new coords in s_c) */
1197 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, c, s_min->s.cg_p, s_c,
1198 constr, top, nrnb, wcycle, -1);
1201 /* Calculate energy for the trial step */
1202 evaluate_energy(fplog, cr,
1203 top_global, s_c, top,
1204 inputrec, nrnb, wcycle, gstat,
1205 vsite, constr, fcd, graph, mdatoms, fr,
1206 mu_tot, enerd, vir, pres, -1, FALSE);
1208 /* Calc derivative along line */
1212 for (i = 0; i < mdatoms->homenr; i++)
1214 for (m = 0; m < DIM; m++)
1216 gpc -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1219 /* Sum the gradient along the line across CPUs */
1222 gmx_sumd(1, &gpc, cr);
1225 /* This is the max amount of increase in energy we tolerate */
1226 tmp = sqrt(GMX_REAL_EPS)*fabs(s_a->epot);
1228 /* Accept the step if the energy is lower, or if it is not significantly higher
1229 * and the line derivative is still negative.
1231 if (s_c->epot < s_a->epot || (gpc < 0 && s_c->epot < (s_a->epot + tmp)))
1234 /* Great, we found a better energy. Increase step for next iteration
1235 * if we are still going down, decrease it otherwise
1239 stepsize *= 1.618034; /* The golden section */
1243 stepsize *= 0.618034; /* 1/golden section */
1248 /* New energy is the same or higher. We will have to do some work
1249 * to find a smaller value in the interval. Take smaller step next time!
1252 stepsize *= 0.618034;
1258 /* OK, if we didn't find a lower value we will have to locate one now - there must
1259 * be one in the interval [a=0,c].
1260 * The same thing is valid here, though: Don't spend dozens of iterations to find
1261 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1262 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1264 * I also have a safeguard for potentially really patological functions so we never
1265 * take more than 20 steps before we give up ...
1267 * If we already found a lower value we just skip this step and continue to the update.
1275 /* Select a new trial point.
1276 * If the derivatives at points a & c have different sign we interpolate to zero,
1277 * otherwise just do a bisection.
1279 if (gpa < 0 && gpc > 0)
1281 b = a + gpa*(a-c)/(gpc-gpa);
1288 /* safeguard if interpolation close to machine accuracy causes errors:
1289 * never go outside the interval
1291 if (b <= a || b >= c)
1296 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
1298 /* Reload the old state */
1299 em_dd_partition_system(fplog, -1, cr, top_global, inputrec,
1300 s_min, top, mdatoms, fr, vsite, constr,
1304 /* Take a trial step to this new point - new coords in s_b */
1305 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC, s_min, b, s_min->s.cg_p, s_b,
1306 constr, top, nrnb, wcycle, -1);
1309 /* Calculate energy for the trial step */
1310 evaluate_energy(fplog, cr,
1311 top_global, s_b, top,
1312 inputrec, nrnb, wcycle, gstat,
1313 vsite, constr, fcd, graph, mdatoms, fr,
1314 mu_tot, enerd, vir, pres, -1, FALSE);
1316 /* p does not change within a step, but since the domain decomposition
1317 * might change, we have to use cg_p of s_b here.
1322 for (i = 0; i < mdatoms->homenr; i++)
1324 for (m = 0; m < DIM; m++)
1326 gpb -= p[i][m]*sf[i][m]; /* f is negative gradient, thus the sign */
1329 /* Sum the gradient along the line across CPUs */
1332 gmx_sumd(1, &gpb, cr);
1337 fprintf(debug, "CGE: EpotA %f EpotB %f EpotC %f gpb %f\n",
1338 s_a->epot, s_b->epot, s_c->epot, gpb);
1341 epot_repl = s_b->epot;
1343 /* Keep one of the intervals based on the value of the derivative at the new point */
1346 /* Replace c endpoint with b */
1347 swap_em_state(s_b, s_c);
1353 /* Replace a endpoint with b */
1354 swap_em_state(s_b, s_a);
1360 * Stop search as soon as we find a value smaller than the endpoints.
1361 * Never run more than 20 steps, no matter what.
1365 while ((epot_repl > s_a->epot || epot_repl > s_c->epot) &&
1368 if (fabs(epot_repl - s_min->epot) < fabs(s_min->epot)*GMX_REAL_EPS ||
1371 /* OK. We couldn't find a significantly lower energy.
1372 * If beta==0 this was steepest descent, and then we give up.
1373 * If not, set beta=0 and restart with steepest descent before quitting.
1383 /* Reset memory before giving up */
1389 /* Select min energy state of A & C, put the best in B.
1391 if (s_c->epot < s_a->epot)
1395 fprintf(debug, "CGE: C (%f) is lower than A (%f), moving C to B\n",
1396 s_c->epot, s_a->epot);
1398 swap_em_state(s_b, s_c);
1406 fprintf(debug, "CGE: A (%f) is lower than C (%f), moving A to B\n",
1407 s_a->epot, s_c->epot);
1409 swap_em_state(s_b, s_a);
1419 fprintf(debug, "CGE: Found a lower energy %f, moving C to B\n",
1422 swap_em_state(s_b, s_c);
1427 /* new search direction */
1428 /* beta = 0 means forget all memory and restart with steepest descents. */
1429 if (nstcg && ((step % nstcg) == 0))
1435 /* s_min->fnorm cannot be zero, because then we would have converged
1439 /* Polak-Ribiere update.
1440 * Change to fnorm2/fnorm2_old for Fletcher-Reeves
1442 beta = pr_beta(cr, &inputrec->opts, mdatoms, top_global, s_min, s_b);
1444 /* Limit beta to prevent oscillations */
1445 if (fabs(beta) > 5.0)
1451 /* update positions */
1452 swap_em_state(s_min, s_b);
1455 /* Print it if necessary */
1460 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
1461 step, s_min->epot, s_min->fnorm/sqrt(state_global->natoms),
1462 s_min->fmax, s_min->a_fmax+1);
1464 /* Store the new (lower) energies */
1465 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1466 mdatoms->tmass, enerd, &s_min->s, inputrec->fepvals, inputrec->expandedvals, s_min->s.box,
1467 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1469 do_log = do_per_step(step, inputrec->nstlog);
1470 do_ene = do_per_step(step, inputrec->nstenergy);
1472 /* Prepare IMD energy record, if bIMD is TRUE. */
1473 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, step, TRUE);
1477 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1479 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
1480 do_log ? fplog : NULL, step, step, eprNORMAL,
1481 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1484 /* Send energies and positions to the IMD client if bIMD is TRUE. */
1485 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
1487 IMD_send_positions(inputrec->imd);
1490 /* Stop when the maximum force lies below tolerance.
1491 * If we have reached machine precision, converged is already set to true.
1493 converged = converged || (s_min->fmax < inputrec->em_tol);
1495 } /* End of the loop */
1497 /* IMD cleanup, if bIMD is TRUE. */
1498 IMD_finalize(inputrec->bIMD, inputrec->imd);
1502 step--; /* we never took that last step in this case */
1505 if (s_min->fmax > inputrec->em_tol)
1509 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
1510 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
1517 /* If we printed energy and/or logfile last step (which was the last step)
1518 * we don't have to do it again, but otherwise print the final values.
1522 /* Write final value to log since we didn't do anything the last step */
1523 print_ebin_header(fplog, step, step, s_min->s.lambda[efptFEP]);
1525 if (!do_ene || !do_log)
1527 /* Write final energy file entries */
1528 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
1529 !do_log ? fplog : NULL, step, step, eprNORMAL,
1530 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1534 /* Print some stuff... */
1537 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
1541 * For accurate normal mode calculation it is imperative that we
1542 * store the last conformation into the full precision binary trajectory.
1544 * However, we should only do it if we did NOT already write this step
1545 * above (which we did if do_x or do_f was true).
1547 do_x = !do_per_step(step, inputrec->nstxout);
1548 do_f = (inputrec->nstfout > 0 && !do_per_step(step, inputrec->nstfout));
1550 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
1551 top_global, inputrec, step,
1552 s_min, state_global, f_global);
1554 fnormn = s_min->fnorm/sqrt(state_global->natoms);
1558 print_converged(stderr, CG, inputrec->em_tol, step, converged, number_steps,
1559 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1560 print_converged(fplog, CG, inputrec->em_tol, step, converged, number_steps,
1561 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
1563 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
1566 finish_em(cr, outf, walltime_accounting, wcycle);
1568 /* To print the actual number of steps we needed somewhere */
1569 walltime_accounting_set_nsteps_done(walltime_accounting, step);
1572 } /* That's all folks */
1575 double do_lbfgs(FILE *fplog, t_commrec *cr,
1576 int nfile, const t_filenm fnm[],
1577 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
1578 int gmx_unused nstglobalcomm,
1579 gmx_vsite_t *vsite, gmx_constr_t constr,
1580 int gmx_unused stepout,
1581 t_inputrec *inputrec,
1582 gmx_mtop_t *top_global, t_fcdata *fcd,
1585 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1586 gmx_edsam_t gmx_unused ed,
1588 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
1589 gmx_membed_t gmx_unused membed,
1590 real gmx_unused cpt_period, real gmx_unused max_hours,
1591 const char gmx_unused *deviceOptions,
1593 unsigned long gmx_unused Flags,
1594 gmx_walltime_accounting_t walltime_accounting)
1596 static const char *LBFGS = "Low-Memory BFGS Minimizer";
1598 gmx_localtop_t *top;
1599 gmx_enerdata_t *enerd;
1601 gmx_global_stat_t gstat;
1604 int ncorr, nmaxcorr, point, cp, neval, nminstep;
1605 double stepsize, gpa, gpb, gpc, tmp, minstep;
1606 real *rho, *alpha, *ff, *xx, *p, *s, *lastx, *lastf, **dx, **dg;
1607 real *xa, *xb, *xc, *fa, *fb, *fc, *xtmp, *ftmp;
1608 real a, b, c, maxdelta, delta;
1609 real diag, Epot0, Epot, EpotA, EpotB, EpotC;
1610 real dgdx, dgdg, sq, yr, beta;
1612 gmx_bool converged, first;
1615 gmx_bool do_log, do_ene, do_x, do_f, foundlower, *frozen;
1617 int start, end, number_steps;
1619 int i, k, m, n, nfmax, gf, step;
1626 gmx_fatal(FARGS, "Cannot do parallel L-BFGS Minimization - yet.\n");
1631 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).");
1634 n = 3*state->natoms;
1635 nmaxcorr = inputrec->nbfgscorr;
1637 /* Allocate memory */
1638 /* Use pointers to real so we dont have to loop over both atoms and
1639 * dimensions all the time...
1640 * x/f are allocated as rvec *, so make new x0/f0 pointers-to-real
1641 * that point to the same memory.
1654 snew(rho, nmaxcorr);
1655 snew(alpha, nmaxcorr);
1658 for (i = 0; i < nmaxcorr; i++)
1664 for (i = 0; i < nmaxcorr; i++)
1673 init_em(fplog, LBFGS, cr, inputrec,
1674 state, top_global, &ems, &top, &f, &f_global,
1675 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
1676 nfile, fnm, &outf, &mdebin, imdport, Flags);
1677 /* Do_lbfgs is not completely updated like do_steep and do_cg,
1678 * so we free some memory again.
1683 xx = (real *)state->x;
1687 end = mdatoms->homenr;
1689 /* Print to log file */
1690 print_em_start(fplog, cr, walltime_accounting, wcycle, LBFGS);
1692 do_log = do_ene = do_x = do_f = TRUE;
1694 /* Max number of steps */
1695 number_steps = inputrec->nsteps;
1697 /* Create a 3*natoms index to tell whether each degree of freedom is frozen */
1699 for (i = start; i < end; i++)
1701 if (mdatoms->cFREEZE)
1703 gf = mdatoms->cFREEZE[i];
1705 for (m = 0; m < DIM; m++)
1707 frozen[3*i+m] = inputrec->opts.nFreeze[gf][m];
1712 sp_header(stderr, LBFGS, inputrec->em_tol, number_steps);
1716 sp_header(fplog, LBFGS, inputrec->em_tol, number_steps);
1721 construct_vsites(vsite, state->x, 1, NULL,
1722 top->idef.iparams, top->idef.il,
1723 fr->ePBC, fr->bMolPBC, cr, state->box);
1726 /* Call the force routine and some auxiliary (neighboursearching etc.) */
1727 /* do_force always puts the charge groups in the box and shifts again
1728 * We do not unshift, so molecules are always whole
1733 evaluate_energy(fplog, cr,
1734 top_global, &ems, top,
1735 inputrec, nrnb, wcycle, gstat,
1736 vsite, constr, fcd, graph, mdatoms, fr,
1737 mu_tot, enerd, vir, pres, -1, TRUE);
1742 /* Copy stuff to the energy bin for easy printing etc. */
1743 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
1744 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
1745 NULL, NULL, vir, pres, NULL, mu_tot, constr);
1747 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
1748 print_ebin(mdoutf_get_fp_ene(outf), TRUE, FALSE, FALSE, fplog, step, step, eprNORMAL,
1749 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
1753 /* This is the starting energy */
1754 Epot = enerd->term[F_EPOT];
1760 /* Set the initial step.
1761 * since it will be multiplied by the non-normalized search direction
1762 * vector (force vector the first time), we scale it by the
1763 * norm of the force.
1768 fprintf(stderr, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1769 fprintf(stderr, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1770 fprintf(stderr, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1771 fprintf(stderr, "\n");
1772 /* and copy to the log file too... */
1773 fprintf(fplog, "Using %d BFGS correction steps.\n\n", nmaxcorr);
1774 fprintf(fplog, " F-max = %12.5e on atom %d\n", fmax, nfmax+1);
1775 fprintf(fplog, " F-Norm = %12.5e\n", fnorm/sqrt(state->natoms));
1776 fprintf(fplog, "\n");
1780 for (i = 0; i < n; i++)
1784 dx[point][i] = ff[i]; /* Initial search direction */
1792 stepsize = 1.0/fnorm;
1794 /* Start the loop over BFGS steps.
1795 * Each successful step is counted, and we continue until
1796 * we either converge or reach the max number of steps.
1801 /* Set the gradient from the force */
1803 for (step = 0; (number_steps < 0 || (number_steps >= 0 && step <= number_steps)) && !converged; step++)
1806 /* Write coordinates if necessary */
1807 do_x = do_per_step(step, inputrec->nstxout);
1808 do_f = do_per_step(step, inputrec->nstfout);
1813 mdof_flags |= MDOF_X;
1818 mdof_flags |= MDOF_F;
1823 mdof_flags |= MDOF_IMD;
1826 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags,
1827 top_global, step, (real)step, state, state, f, f);
1829 /* Do the linesearching in the direction dx[point][0..(n-1)] */
1831 /* pointer to current direction - point=0 first time here */
1834 /* calculate line gradient */
1835 for (gpa = 0, i = 0; i < n; i++)
1840 /* Calculate minimum allowed stepsize, before the average (norm)
1841 * relative change in coordinate is smaller than precision
1843 for (minstep = 0, i = 0; i < n; i++)
1853 minstep = GMX_REAL_EPS/sqrt(minstep/n);
1855 if (stepsize < minstep)
1861 /* Store old forces and coordinates */
1862 for (i = 0; i < n; i++)
1871 for (i = 0; i < n; i++)
1876 /* Take a step downhill.
1877 * In theory, we should minimize the function along this direction.
1878 * That is quite possible, but it turns out to take 5-10 function evaluations
1879 * for each line. However, we dont really need to find the exact minimum -
1880 * it is much better to start a new BFGS step in a modified direction as soon
1881 * as we are close to it. This will save a lot of energy evaluations.
1883 * In practice, we just try to take a single step.
1884 * If it worked (i.e. lowered the energy), we increase the stepsize but
1885 * the continue straight to the next BFGS step without trying to find any minimum.
1886 * If it didn't work (higher energy), there must be a minimum somewhere between
1887 * the old position and the new one.
1889 * Due to the finite numerical accuracy, it turns out that it is a good idea
1890 * to even accept a SMALL increase in energy, if the derivative is still downhill.
1891 * This leads to lower final energies in the tests I've done. / Erik
1896 c = a + stepsize; /* reference position along line is zero */
1898 /* Check stepsize first. We do not allow displacements
1899 * larger than emstep.
1905 for (i = 0; i < n; i++)
1908 if (delta > maxdelta)
1913 if (maxdelta > inputrec->em_stepsize)
1918 while (maxdelta > inputrec->em_stepsize);
1920 /* Take a trial step */
1921 for (i = 0; i < n; i++)
1923 xc[i] = lastx[i] + c*s[i];
1927 /* Calculate energy for the trial step */
1928 ems.s.x = (rvec *)xc;
1930 evaluate_energy(fplog, cr,
1931 top_global, &ems, top,
1932 inputrec, nrnb, wcycle, gstat,
1933 vsite, constr, fcd, graph, mdatoms, fr,
1934 mu_tot, enerd, vir, pres, step, FALSE);
1937 /* Calc derivative along line */
1938 for (gpc = 0, i = 0; i < n; i++)
1940 gpc -= s[i]*fc[i]; /* f is negative gradient, thus the sign */
1942 /* Sum the gradient along the line across CPUs */
1945 gmx_sumd(1, &gpc, cr);
1948 /* This is the max amount of increase in energy we tolerate */
1949 tmp = sqrt(GMX_REAL_EPS)*fabs(EpotA);
1951 /* Accept the step if the energy is lower, or if it is not significantly higher
1952 * and the line derivative is still negative.
1954 if (EpotC < EpotA || (gpc < 0 && EpotC < (EpotA+tmp)))
1957 /* Great, we found a better energy. Increase step for next iteration
1958 * if we are still going down, decrease it otherwise
1962 stepsize *= 1.618034; /* The golden section */
1966 stepsize *= 0.618034; /* 1/golden section */
1971 /* New energy is the same or higher. We will have to do some work
1972 * to find a smaller value in the interval. Take smaller step next time!
1975 stepsize *= 0.618034;
1978 /* OK, if we didn't find a lower value we will have to locate one now - there must
1979 * be one in the interval [a=0,c].
1980 * The same thing is valid here, though: Don't spend dozens of iterations to find
1981 * the line minimum. We try to interpolate based on the derivative at the endpoints,
1982 * and only continue until we find a lower value. In most cases this means 1-2 iterations.
1984 * I also have a safeguard for potentially really patological functions so we never
1985 * take more than 20 steps before we give up ...
1987 * If we already found a lower value we just skip this step and continue to the update.
1996 /* Select a new trial point.
1997 * If the derivatives at points a & c have different sign we interpolate to zero,
1998 * otherwise just do a bisection.
2001 if (gpa < 0 && gpc > 0)
2003 b = a + gpa*(a-c)/(gpc-gpa);
2010 /* safeguard if interpolation close to machine accuracy causes errors:
2011 * never go outside the interval
2013 if (b <= a || b >= c)
2018 /* Take a trial step */
2019 for (i = 0; i < n; i++)
2021 xb[i] = lastx[i] + b*s[i];
2025 /* Calculate energy for the trial step */
2026 ems.s.x = (rvec *)xb;
2028 evaluate_energy(fplog, cr,
2029 top_global, &ems, top,
2030 inputrec, nrnb, wcycle, gstat,
2031 vsite, constr, fcd, graph, mdatoms, fr,
2032 mu_tot, enerd, vir, pres, step, FALSE);
2037 for (gpb = 0, i = 0; i < n; i++)
2039 gpb -= s[i]*fb[i]; /* f is negative gradient, thus the sign */
2042 /* Sum the gradient along the line across CPUs */
2045 gmx_sumd(1, &gpb, cr);
2048 /* Keep one of the intervals based on the value of the derivative at the new point */
2051 /* Replace c endpoint with b */
2055 /* swap coord pointers b/c */
2065 /* Replace a endpoint with b */
2069 /* swap coord pointers a/b */
2079 * Stop search as soon as we find a value smaller than the endpoints,
2080 * or if the tolerance is below machine precision.
2081 * Never run more than 20 steps, no matter what.
2085 while ((EpotB > EpotA || EpotB > EpotC) && (nminstep < 20));
2087 if (fabs(EpotB-Epot0) < GMX_REAL_EPS || nminstep >= 20)
2089 /* OK. We couldn't find a significantly lower energy.
2090 * If ncorr==0 this was steepest descent, and then we give up.
2091 * If not, reset memory to restart as steepest descent before quitting.
2103 /* Search in gradient direction */
2104 for (i = 0; i < n; i++)
2106 dx[point][i] = ff[i];
2108 /* Reset stepsize */
2109 stepsize = 1.0/fnorm;
2114 /* Select min energy state of A & C, put the best in xx/ff/Epot
2120 for (i = 0; i < n; i++)
2131 for (i = 0; i < n; i++)
2145 for (i = 0; i < n; i++)
2153 /* Update the memory information, and calculate a new
2154 * approximation of the inverse hessian
2157 /* Have new data in Epot, xx, ff */
2158 if (ncorr < nmaxcorr)
2163 for (i = 0; i < n; i++)
2165 dg[point][i] = lastf[i]-ff[i];
2166 dx[point][i] *= stepsize;
2171 for (i = 0; i < n; i++)
2173 dgdg += dg[point][i]*dg[point][i];
2174 dgdx += dg[point][i]*dx[point][i];
2179 rho[point] = 1.0/dgdx;
2182 if (point >= nmaxcorr)
2188 for (i = 0; i < n; i++)
2195 /* Recursive update. First go back over the memory points */
2196 for (k = 0; k < ncorr; k++)
2205 for (i = 0; i < n; i++)
2207 sq += dx[cp][i]*p[i];
2210 alpha[cp] = rho[cp]*sq;
2212 for (i = 0; i < n; i++)
2214 p[i] -= alpha[cp]*dg[cp][i];
2218 for (i = 0; i < n; i++)
2223 /* And then go forward again */
2224 for (k = 0; k < ncorr; k++)
2227 for (i = 0; i < n; i++)
2229 yr += p[i]*dg[cp][i];
2233 beta = alpha[cp]-beta;
2235 for (i = 0; i < n; i++)
2237 p[i] += beta*dx[cp][i];
2247 for (i = 0; i < n; i++)
2251 dx[point][i] = p[i];
2261 /* Test whether the convergence criterion is met */
2262 get_f_norm_max(cr, &(inputrec->opts), mdatoms, f, &fnorm, &fmax, &nfmax);
2264 /* Print it if necessary */
2269 fprintf(stderr, "\rStep %d, Epot=%12.6e, Fnorm=%9.3e, Fmax=%9.3e (atom %d)\n",
2270 step, Epot, fnorm/sqrt(state->natoms), fmax, nfmax+1);
2272 /* Store the new (lower) energies */
2273 upd_mdebin(mdebin, FALSE, FALSE, (double)step,
2274 mdatoms->tmass, enerd, state, inputrec->fepvals, inputrec->expandedvals, state->box,
2275 NULL, NULL, vir, pres, NULL, mu_tot, constr);
2276 do_log = do_per_step(step, inputrec->nstlog);
2277 do_ene = do_per_step(step, inputrec->nstenergy);
2280 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2282 print_ebin(mdoutf_get_fp_ene(outf), do_ene, FALSE, FALSE,
2283 do_log ? fplog : NULL, step, step, eprNORMAL,
2284 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2287 /* Send x and E to IMD client, if bIMD is TRUE. */
2288 if (do_IMD(inputrec->bIMD, step, cr, TRUE, state->box, state->x, inputrec, 0, wcycle) && MASTER(cr))
2290 IMD_send_positions(inputrec->imd);
2293 /* Stop when the maximum force lies below tolerance.
2294 * If we have reached machine precision, converged is already set to true.
2297 converged = converged || (fmax < inputrec->em_tol);
2299 } /* End of the loop */
2301 /* IMD cleanup, if bIMD is TRUE. */
2302 IMD_finalize(inputrec->bIMD, inputrec->imd);
2306 step--; /* we never took that last step in this case */
2309 if (fmax > inputrec->em_tol)
2313 warn_step(stderr, inputrec->em_tol, step-1 == number_steps, FALSE);
2314 warn_step(fplog, inputrec->em_tol, step-1 == number_steps, FALSE);
2319 /* If we printed energy and/or logfile last step (which was the last step)
2320 * we don't have to do it again, but otherwise print the final values.
2322 if (!do_log) /* Write final value to log since we didn't do anythin last step */
2324 print_ebin_header(fplog, step, step, state->lambda[efptFEP]);
2326 if (!do_ene || !do_log) /* Write final energy file entries */
2328 print_ebin(mdoutf_get_fp_ene(outf), !do_ene, FALSE, FALSE,
2329 !do_log ? fplog : NULL, step, step, eprNORMAL,
2330 TRUE, mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2333 /* Print some stuff... */
2336 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2340 * For accurate normal mode calculation it is imperative that we
2341 * store the last conformation into the full precision binary trajectory.
2343 * However, we should only do it if we did NOT already write this step
2344 * above (which we did if do_x or do_f was true).
2346 do_x = !do_per_step(step, inputrec->nstxout);
2347 do_f = !do_per_step(step, inputrec->nstfout);
2348 write_em_traj(fplog, cr, outf, do_x, do_f, ftp2fn(efSTO, nfile, fnm),
2349 top_global, inputrec, step,
2354 print_converged(stderr, LBFGS, inputrec->em_tol, step, converged,
2355 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2356 print_converged(fplog, LBFGS, inputrec->em_tol, step, converged,
2357 number_steps, Epot, fmax, nfmax, fnorm/sqrt(state->natoms));
2359 fprintf(fplog, "\nPerformed %d energy evaluations in total.\n", neval);
2362 finish_em(cr, outf, walltime_accounting, wcycle);
2364 /* To print the actual number of steps we needed somewhere */
2365 walltime_accounting_set_nsteps_done(walltime_accounting, step);
2368 } /* That's all folks */
2371 double do_steep(FILE *fplog, t_commrec *cr,
2372 int nfile, const t_filenm fnm[],
2373 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2374 int gmx_unused nstglobalcomm,
2375 gmx_vsite_t *vsite, gmx_constr_t constr,
2376 int gmx_unused stepout,
2377 t_inputrec *inputrec,
2378 gmx_mtop_t *top_global, t_fcdata *fcd,
2379 t_state *state_global,
2381 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2382 gmx_edsam_t gmx_unused ed,
2384 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2385 gmx_membed_t gmx_unused membed,
2386 real gmx_unused cpt_period, real gmx_unused max_hours,
2387 const char gmx_unused *deviceOptions,
2389 unsigned long gmx_unused Flags,
2390 gmx_walltime_accounting_t walltime_accounting)
2392 const char *SD = "Steepest Descents";
2393 em_state_t *s_min, *s_try;
2395 gmx_localtop_t *top;
2396 gmx_enerdata_t *enerd;
2398 gmx_global_stat_t gstat;
2400 real stepsize, constepsize;
2404 gmx_bool bDone, bAbort, do_x, do_f;
2409 int steps_accepted = 0;
2413 s_min = init_em_state();
2414 s_try = init_em_state();
2416 /* Init em and store the local state in s_try */
2417 init_em(fplog, SD, cr, inputrec,
2418 state_global, top_global, s_try, &top, &f, &f_global,
2419 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2420 nfile, fnm, &outf, &mdebin, imdport, Flags);
2422 /* Print to log file */
2423 print_em_start(fplog, cr, walltime_accounting, wcycle, SD);
2425 /* Set variables for stepsize (in nm). This is the largest
2426 * step that we are going to make in any direction.
2428 ustep = inputrec->em_stepsize;
2431 /* Max number of steps */
2432 nsteps = inputrec->nsteps;
2436 /* Print to the screen */
2437 sp_header(stderr, SD, inputrec->em_tol, nsteps);
2441 sp_header(fplog, SD, inputrec->em_tol, nsteps);
2444 /**** HERE STARTS THE LOOP ****
2445 * count is the counter for the number of steps
2446 * bDone will be TRUE when the minimization has converged
2447 * bAbort will be TRUE when nsteps steps have been performed or when
2448 * the stepsize becomes smaller than is reasonable for machine precision
2453 while (!bDone && !bAbort)
2455 bAbort = (nsteps >= 0) && (count == nsteps);
2457 /* set new coordinates, except for first step */
2460 do_em_step(cr, inputrec, mdatoms, fr->bMolPBC,
2461 s_min, stepsize, s_min->f, s_try,
2462 constr, top, nrnb, wcycle, count);
2465 evaluate_energy(fplog, cr,
2466 top_global, s_try, top,
2467 inputrec, nrnb, wcycle, gstat,
2468 vsite, constr, fcd, graph, mdatoms, fr,
2469 mu_tot, enerd, vir, pres, count, count == 0);
2473 print_ebin_header(fplog, count, count, s_try->s.lambda[efptFEP]);
2478 s_min->epot = s_try->epot + 1;
2481 /* Print it if necessary */
2486 fprintf(stderr, "Step=%5d, Dmax= %6.1e nm, Epot= %12.5e Fmax= %11.5e, atom= %d%c",
2487 count, ustep, s_try->epot, s_try->fmax, s_try->a_fmax+1,
2488 (s_try->epot < s_min->epot) ? '\n' : '\r');
2491 if (s_try->epot < s_min->epot)
2493 /* Store the new (lower) energies */
2494 upd_mdebin(mdebin, FALSE, FALSE, (double)count,
2495 mdatoms->tmass, enerd, &s_try->s, inputrec->fepvals, inputrec->expandedvals,
2496 s_try->s.box, NULL, NULL, vir, pres, NULL, mu_tot, constr);
2498 /* Prepare IMD energy record, if bIMD is TRUE. */
2499 IMD_fill_energy_record(inputrec->bIMD, inputrec->imd, enerd, count, TRUE);
2501 print_ebin(mdoutf_get_fp_ene(outf), TRUE,
2502 do_per_step(steps_accepted, inputrec->nstdisreout),
2503 do_per_step(steps_accepted, inputrec->nstorireout),
2504 fplog, count, count, eprNORMAL, TRUE,
2505 mdebin, fcd, &(top_global->groups), &(inputrec->opts));
2510 /* Now if the new energy is smaller than the previous...
2511 * or if this is the first step!
2512 * or if we did random steps!
2515 if ( (count == 0) || (s_try->epot < s_min->epot) )
2519 /* Test whether the convergence criterion is met... */
2520 bDone = (s_try->fmax < inputrec->em_tol);
2522 /* Copy the arrays for force, positions and energy */
2523 /* The 'Min' array always holds the coords and forces of the minimal
2525 swap_em_state(s_min, s_try);
2531 /* Write to trn, if necessary */
2532 do_x = do_per_step(steps_accepted, inputrec->nstxout);
2533 do_f = do_per_step(steps_accepted, inputrec->nstfout);
2534 write_em_traj(fplog, cr, outf, do_x, do_f, NULL,
2535 top_global, inputrec, count,
2536 s_min, state_global, f_global);
2540 /* If energy is not smaller make the step smaller... */
2543 if (DOMAINDECOMP(cr) && s_min->s.ddp_count != cr->dd->ddp_count)
2545 /* Reload the old state */
2546 em_dd_partition_system(fplog, count, cr, top_global, inputrec,
2547 s_min, top, mdatoms, fr, vsite, constr,
2552 /* Determine new step */
2553 stepsize = ustep/s_min->fmax;
2555 /* Check if stepsize is too small, with 1 nm as a characteristic length */
2557 if (count == nsteps || ustep < 1e-12)
2559 if (count == nsteps || ustep < 1e-6)
2564 warn_step(stderr, inputrec->em_tol, count == nsteps, constr != NULL);
2565 warn_step(fplog, inputrec->em_tol, count == nsteps, constr != NULL);
2570 /* Send IMD energies and positions, if bIMD is TRUE. */
2571 if (do_IMD(inputrec->bIMD, count, cr, TRUE, state_global->box, state_global->x, inputrec, 0, wcycle) && MASTER(cr))
2573 IMD_send_positions(inputrec->imd);
2577 } /* End of the loop */
2579 /* IMD cleanup, if bIMD is TRUE. */
2580 IMD_finalize(inputrec->bIMD, inputrec->imd);
2582 /* Print some data... */
2585 fprintf(stderr, "\nwriting lowest energy coordinates.\n");
2587 write_em_traj(fplog, cr, outf, TRUE, inputrec->nstfout, ftp2fn(efSTO, nfile, fnm),
2588 top_global, inputrec, count,
2589 s_min, state_global, f_global);
2591 fnormn = s_min->fnorm/sqrt(state_global->natoms);
2595 print_converged(stderr, SD, inputrec->em_tol, count, bDone, nsteps,
2596 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2597 print_converged(fplog, SD, inputrec->em_tol, count, bDone, nsteps,
2598 s_min->epot, s_min->fmax, s_min->a_fmax, fnormn);
2601 finish_em(cr, outf, walltime_accounting, wcycle);
2603 /* To print the actual number of steps we needed somewhere */
2604 inputrec->nsteps = count;
2606 walltime_accounting_set_nsteps_done(walltime_accounting, count);
2609 } /* That's all folks */
2612 double do_nm(FILE *fplog, t_commrec *cr,
2613 int nfile, const t_filenm fnm[],
2614 const output_env_t gmx_unused oenv, gmx_bool bVerbose, gmx_bool gmx_unused bCompact,
2615 int gmx_unused nstglobalcomm,
2616 gmx_vsite_t *vsite, gmx_constr_t constr,
2617 int gmx_unused stepout,
2618 t_inputrec *inputrec,
2619 gmx_mtop_t *top_global, t_fcdata *fcd,
2620 t_state *state_global,
2622 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
2623 gmx_edsam_t gmx_unused ed,
2625 int gmx_unused repl_ex_nst, int gmx_unused repl_ex_nex, int gmx_unused repl_ex_seed,
2626 gmx_membed_t gmx_unused membed,
2627 real gmx_unused cpt_period, real gmx_unused max_hours,
2628 const char gmx_unused *deviceOptions,
2630 unsigned long gmx_unused Flags,
2631 gmx_walltime_accounting_t walltime_accounting)
2633 const char *NM = "Normal Mode Analysis";
2635 int natoms, atom, d;
2638 gmx_localtop_t *top;
2639 gmx_enerdata_t *enerd;
2641 gmx_global_stat_t gstat;
2643 real t, t0, lambda, lam0;
2648 gmx_bool bSparse; /* use sparse matrix storage format */
2650 gmx_sparsematrix_t * sparse_matrix = NULL;
2651 real * full_matrix = NULL;
2652 em_state_t * state_work;
2654 /* added with respect to mdrun */
2655 int i, j, k, row, col;
2656 real der_range = 10.0*sqrt(GMX_REAL_EPS);
2662 gmx_fatal(FARGS, "Constraints present with Normal Mode Analysis, this combination is not supported");
2665 state_work = init_em_state();
2667 /* Init em and store the local state in state_minimum */
2668 init_em(fplog, NM, cr, inputrec,
2669 state_global, top_global, state_work, &top,
2671 nrnb, mu_tot, fr, &enerd, &graph, mdatoms, &gstat, vsite, constr,
2672 nfile, fnm, &outf, NULL, imdport, Flags);
2674 natoms = top_global->natoms;
2682 "NOTE: This version of Gromacs has been compiled in single precision,\n"
2683 " which MIGHT not be accurate enough for normal mode analysis.\n"
2684 " Gromacs now uses sparse matrix storage, so the memory requirements\n"
2685 " are fairly modest even if you recompile in double precision.\n\n");
2689 /* Check if we can/should use sparse storage format.
2691 * Sparse format is only useful when the Hessian itself is sparse, which it
2692 * will be when we use a cutoff.
2693 * For small systems (n<1000) it is easier to always use full matrix format, though.
2695 if (EEL_FULL(fr->eeltype) || fr->rlist == 0.0)
2697 md_print_info(cr, fplog, "Non-cutoff electrostatics used, forcing full Hessian format.\n");
2700 else if (top_global->natoms < 1000)
2702 md_print_info(cr, fplog, "Small system size (N=%d), using full Hessian format.\n", top_global->natoms);
2707 md_print_info(cr, fplog, "Using compressed symmetric sparse Hessian format.\n");
2713 sz = DIM*top_global->natoms;
2715 fprintf(stderr, "Allocating Hessian memory...\n\n");
2719 sparse_matrix = gmx_sparsematrix_init(sz);
2720 sparse_matrix->compressed_symmetric = TRUE;
2724 snew(full_matrix, sz*sz);
2728 /* Initial values */
2729 t0 = inputrec->init_t;
2730 lam0 = inputrec->fepvals->init_lambda;
2738 /* Write start time and temperature */
2739 print_em_start(fplog, cr, walltime_accounting, wcycle, NM);
2741 /* fudge nr of steps to nr of atoms */
2742 inputrec->nsteps = natoms*2;
2746 fprintf(stderr, "starting normal mode calculation '%s'\n%d steps.\n\n",
2747 *(top_global->name), (int)inputrec->nsteps);
2750 nnodes = cr->nnodes;
2752 /* Make evaluate_energy do a single node force calculation */
2754 evaluate_energy(fplog, cr,
2755 top_global, state_work, top,
2756 inputrec, nrnb, wcycle, gstat,
2757 vsite, constr, fcd, graph, mdatoms, fr,
2758 mu_tot, enerd, vir, pres, -1, TRUE);
2759 cr->nnodes = nnodes;
2761 /* if forces are not small, warn user */
2762 get_state_f_norm_max(cr, &(inputrec->opts), mdatoms, state_work);
2764 md_print_info(cr, fplog, "Maximum force:%12.5e\n", state_work->fmax);
2765 if (state_work->fmax > 1.0e-3)
2767 md_print_info(cr, fplog,
2768 "The force is probably not small enough to "
2769 "ensure that you are at a minimum.\n"
2770 "Be aware that negative eigenvalues may occur\n"
2771 "when the resulting matrix is diagonalized.\n\n");
2774 /***********************************************************
2776 * Loop over all pairs in matrix
2778 * do_force called twice. Once with positive and
2779 * once with negative displacement
2781 ************************************************************/
2783 /* Steps are divided one by one over the nodes */
2784 for (atom = cr->nodeid; atom < natoms; atom += nnodes)
2787 for (d = 0; d < DIM; d++)
2789 x_min = state_work->s.x[atom][d];
2791 state_work->s.x[atom][d] = x_min - der_range;
2793 /* Make evaluate_energy do a single node force calculation */
2795 evaluate_energy(fplog, cr,
2796 top_global, state_work, top,
2797 inputrec, nrnb, wcycle, gstat,
2798 vsite, constr, fcd, graph, mdatoms, fr,
2799 mu_tot, enerd, vir, pres, atom*2, FALSE);
2801 for (i = 0; i < natoms; i++)
2803 copy_rvec(state_work->f[i], fneg[i]);
2806 state_work->s.x[atom][d] = x_min + der_range;
2808 evaluate_energy(fplog, cr,
2809 top_global, state_work, top,
2810 inputrec, nrnb, wcycle, gstat,
2811 vsite, constr, fcd, graph, mdatoms, fr,
2812 mu_tot, enerd, vir, pres, atom*2+1, FALSE);
2813 cr->nnodes = nnodes;
2815 /* x is restored to original */
2816 state_work->s.x[atom][d] = x_min;
2818 for (j = 0; j < natoms; j++)
2820 for (k = 0; (k < DIM); k++)
2823 -(state_work->f[j][k] - fneg[j][k])/(2*der_range);
2831 #define mpi_type MPI_DOUBLE
2833 #define mpi_type MPI_FLOAT
2835 MPI_Send(dfdx[0], natoms*DIM, mpi_type, MASTERNODE(cr), cr->nodeid,
2836 cr->mpi_comm_mygroup);
2841 for (node = 0; (node < nnodes && atom+node < natoms); node++)
2847 MPI_Recv(dfdx[0], natoms*DIM, mpi_type, node, node,
2848 cr->mpi_comm_mygroup, &stat);
2853 row = (atom + node)*DIM + d;
2855 for (j = 0; j < natoms; j++)
2857 for (k = 0; k < DIM; k++)
2863 if (col >= row && dfdx[j][k] != 0.0)
2865 gmx_sparsematrix_increment_value(sparse_matrix,
2866 row, col, dfdx[j][k]);
2871 full_matrix[row*sz+col] = dfdx[j][k];
2878 if (bVerbose && fplog)
2883 /* write progress */
2884 if (MASTER(cr) && bVerbose)
2886 fprintf(stderr, "\rFinished step %d out of %d",
2887 min(atom+nnodes, natoms), natoms);
2894 fprintf(stderr, "\n\nWriting Hessian...\n");
2895 gmx_mtxio_write(ftp2fn(efMTX, nfile, fnm), sz, sz, full_matrix, sparse_matrix);
2898 finish_em(cr, outf, walltime_accounting, wcycle);
2900 walltime_accounting_set_nsteps_done(walltime_accounting, natoms*2);